Actual source code: mpiov.c
petsc-3.13.6 2020-09-29
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;
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,owner;
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;
370: PetscMPIInt rank,owner;
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 received (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,proc = 0;
471: PetscInt M,i,j,k,**rbuf,row,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 received 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 received 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;
1271: if (ismax != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This routine only works when all processes have ismax=1");
1273: PetscObjectGetComm((PetscObject)C,&comm);
1274: size = c->size;
1275: rank = c->rank;
1277: ISSorted(iscol[0],&iscolsorted);
1278: ISSorted(isrow[0],&isrowsorted);
1279: ISGetIndices(isrow[0],&irow);
1280: ISGetLocalSize(isrow[0],&nrow);
1281: if (allcolumns) {
1282: icol = NULL;
1283: ncol = C->cmap->N;
1284: } else {
1285: ISGetIndices(iscol[0],&icol);
1286: ISGetLocalSize(iscol[0],&ncol);
1287: }
1289: if (scall == MAT_INITIAL_MATRIX) {
1290: PetscInt *sbuf2_i,*cworkA,lwrite,ctmp;
1292: /* Get some new tags to keep the communication clean */
1293: tag1 = ((PetscObject)C)->tag;
1294: PetscObjectGetNewTag((PetscObject)C,&tag2);
1295: PetscObjectGetNewTag((PetscObject)C,&tag3);
1297: /* evaluate communication - mesg to who, length of mesg, and buffer space
1298: required. Based on this, buffers are allocated, and data copied into them */
1299: PetscCalloc2(size,&w1,size,&w2);
1300: PetscMalloc1(nrow,&row2proc);
1302: /* w1[proc] = num of rows owned by proc -- to be requested */
1303: proc = 0;
1304: nrqs = 0; /* num of outgoing messages */
1305: for (j=0; j<nrow; j++) {
1306: row = irow[j];
1307: if (!isrowsorted) proc = 0;
1308: while (row >= C->rmap->range[proc+1]) proc++;
1309: w1[proc]++;
1310: row2proc[j] = proc; /* map row index to proc */
1312: if (proc != rank && !w2[proc]) {
1313: w2[proc] = 1; nrqs++;
1314: }
1315: }
1316: w1[rank] = 0; /* rows owned by self will not be requested */
1318: PetscMalloc1(nrqs+1,&pa); /*(proc -array)*/
1319: for (proc=0,j=0; proc<size; proc++) {
1320: if (w1[proc]) { pa[j++] = proc;}
1321: }
1323: /* Each message would have a header = 1 + 2*(num of IS) + data (here,num of IS = 1) */
1324: msz = 0; /* total mesg length (for all procs) */
1325: for (i=0; i<nrqs; i++) {
1326: proc = pa[i];
1327: w1[proc] += 3;
1328: msz += w1[proc];
1329: }
1330: PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);
1332: /* Determine nrqr, the number of messages to expect, their lengths, from from-ids */
1333: /* if w2[proc]=1, a message of length w1[proc] will be sent to proc; */
1334: PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);
1336: /* Input: nrqs: nsend; nrqr: nrecv; w1: msg length to be sent;
1337: Output: onodes1: recv node-ids; olengths1: corresponding recv message length */
1338: PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);
1340: /* Now post the Irecvs corresponding to these messages */
1341: PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);
1343: PetscFree(onodes1);
1344: PetscFree(olengths1);
1346: /* Allocate Memory for outgoing messages */
1347: PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);
1348: PetscArrayzero(sbuf1,size);
1349: PetscArrayzero(ptr,size);
1351: /* subf1[pa[0]] = tmp, subf1[pa[i]] = subf1[pa[i-1]] + w1[pa[i-1]] */
1352: iptr = tmp;
1353: for (i=0; i<nrqs; i++) {
1354: proc = pa[i];
1355: sbuf1[proc] = iptr;
1356: iptr += w1[proc];
1357: }
1359: /* Form the outgoing messages */
1360: /* Initialize the header space */
1361: for (i=0; i<nrqs; i++) {
1362: proc = pa[i];
1363: PetscArrayzero(sbuf1[proc],3);
1364: ptr[proc] = sbuf1[proc] + 3;
1365: }
1367: /* Parse the isrow and copy data into outbuf */
1368: PetscArrayzero(ctr,size);
1369: for (j=0; j<nrow; j++) { /* parse the indices of each IS */
1370: proc = row2proc[j];
1371: if (proc != rank) { /* copy to the outgoing buf*/
1372: *ptr[proc] = irow[j];
1373: ctr[proc]++; ptr[proc]++;
1374: }
1375: }
1377: /* Update the headers for the current IS */
1378: for (j=0; j<size; j++) { /* Can Optimise this loop too */
1379: if ((ctr_j = ctr[j])) {
1380: sbuf1_j = sbuf1[j];
1381: k = ++sbuf1_j[0];
1382: sbuf1_j[2*k] = ctr_j;
1383: sbuf1_j[2*k-1] = 0;
1384: }
1385: }
1387: /* Now post the sends */
1388: PetscMalloc1(nrqs+1,&s_waits1);
1389: for (i=0; i<nrqs; ++i) {
1390: proc = pa[i];
1391: MPI_Isend(sbuf1[proc],w1[proc],MPIU_INT,proc,tag1,comm,s_waits1+i);
1392: }
1394: /* Post Receives to capture the buffer size */
1395: PetscMalloc4(nrqs+1,&r_status2,nrqr+1,&s_waits2,nrqs+1,&r_waits2,nrqr+1,&s_status2);
1396: PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);
1398: rbuf2[0] = tmp + msz;
1399: for (i=1; i<nrqs; ++i) rbuf2[i] = rbuf2[i-1] + w1[pa[i-1]];
1401: for (i=0; i<nrqs; ++i) {
1402: proc = pa[i];
1403: MPI_Irecv(rbuf2[i],w1[proc],MPIU_INT,proc,tag2,comm,r_waits2+i);
1404: }
1406: PetscFree2(w1,w2);
1408: /* Send to other procs the buf size they should allocate */
1409: /* Receive messages*/
1410: PetscMalloc1(nrqr+1,&r_status1);
1411: PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);
1413: MPI_Waitall(nrqr,r_waits1,r_status1);
1414: for (i=0; i<nrqr; ++i) {
1415: req_size[i] = 0;
1416: rbuf1_i = rbuf1[i];
1417: start = 2*rbuf1_i[0] + 1;
1418: MPI_Get_count(r_status1+i,MPIU_INT,&end);
1419: PetscMalloc1(end+1,&sbuf2[i]);
1420: sbuf2_i = sbuf2[i];
1421: for (j=start; j<end; j++) {
1422: k = rbuf1_i[j] - rstart;
1423: ncols = ai[k+1] - ai[k] + bi[k+1] - bi[k];
1424: sbuf2_i[j] = ncols;
1425: req_size[i] += ncols;
1426: }
1427: req_source1[i] = r_status1[i].MPI_SOURCE;
1429: /* form the header */
1430: sbuf2_i[0] = req_size[i];
1431: for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];
1433: MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);
1434: }
1436: PetscFree(r_status1);
1437: PetscFree(r_waits1);
1439: /* rbuf2 is received, Post recv column indices a->j */
1440: MPI_Waitall(nrqs,r_waits2,r_status2);
1442: PetscMalloc4(nrqs+1,&r_waits3,nrqr+1,&s_waits3,nrqs+1,&r_status3,nrqr+1,&s_status3);
1443: for (i=0; i<nrqs; ++i) {
1444: PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);
1445: req_source2[i] = r_status2[i].MPI_SOURCE;
1446: MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);
1447: }
1449: /* Wait on sends1 and sends2 */
1450: PetscMalloc1(nrqs+1,&s_status1);
1451: MPI_Waitall(nrqs,s_waits1,s_status1);
1452: PetscFree(s_waits1);
1453: PetscFree(s_status1);
1455: MPI_Waitall(nrqr,s_waits2,s_status2);
1456: PetscFree4(r_status2,s_waits2,r_waits2,s_status2);
1458: /* Now allocate sending buffers for a->j, and send them off */
1459: PetscMalloc1(nrqr+1,&sbuf_aj);
1460: for (i=0,j=0; i<nrqr; i++) j += req_size[i];
1461: PetscMalloc1(j+1,&sbuf_aj[0]);
1462: for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
1464: for (i=0; i<nrqr; i++) { /* for each requested message */
1465: rbuf1_i = rbuf1[i];
1466: sbuf_aj_i = sbuf_aj[i];
1467: ct1 = 2*rbuf1_i[0] + 1;
1468: ct2 = 0;
1469: /* max1=rbuf1_i[0]; if (max1 != 1) SETERRQ1(PETSC_COMM_SELF,0,"max1 %d != 1",max1); */
1471: kmax = rbuf1[i][2];
1472: for (k=0; k<kmax; k++,ct1++) { /* for each row */
1473: row = rbuf1_i[ct1] - rstart;
1474: nzA = ai[row+1] - ai[row];
1475: nzB = bi[row+1] - bi[row];
1476: ncols = nzA + nzB;
1477: cworkA = aj + ai[row]; cworkB = bj + bi[row];
1479: /* load the column indices for this row into cols*/
1480: cols = sbuf_aj_i + ct2;
1482: lwrite = 0;
1483: for (l=0; l<nzB; l++) {
1484: if ((ctmp = bmap[cworkB[l]]) < cstart) cols[lwrite++] = ctmp;
1485: }
1486: for (l=0; l<nzA; l++) cols[lwrite++] = cstart + cworkA[l];
1487: for (l=0; l<nzB; l++) {
1488: if ((ctmp = bmap[cworkB[l]]) >= cend) cols[lwrite++] = ctmp;
1489: }
1491: ct2 += ncols;
1492: }
1493: MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);
1494: }
1496: /* create column map (cmap): global col of C -> local col of submat */
1497: #if defined(PETSC_USE_CTABLE)
1498: if (!allcolumns) {
1499: PetscTableCreate(ncol+1,C->cmap->N+1,&cmap);
1500: PetscCalloc1(C->cmap->n,&cmap_loc);
1501: for (j=0; j<ncol; j++) { /* use array cmap_loc[] for local col indices */
1502: if (icol[j] >= cstart && icol[j] <cend) {
1503: cmap_loc[icol[j] - cstart] = j+1;
1504: } else { /* use PetscTable for non-local col indices */
1505: PetscTableAdd(cmap,icol[j]+1,j+1,INSERT_VALUES);
1506: }
1507: }
1508: } else {
1509: cmap = NULL;
1510: cmap_loc = NULL;
1511: }
1512: PetscCalloc1(C->rmap->n,&rmap_loc);
1513: #else
1514: if (!allcolumns) {
1515: PetscCalloc1(C->cmap->N,&cmap);
1516: for (j=0; j<ncol; j++) cmap[icol[j]] = j+1;
1517: } else {
1518: cmap = NULL;
1519: }
1520: #endif
1522: /* Create lens for MatSeqAIJSetPreallocation() */
1523: PetscCalloc1(nrow,&lens);
1525: /* Compute lens from local part of C */
1526: for (j=0; j<nrow; j++) {
1527: row = irow[j];
1528: proc = row2proc[j];
1529: if (proc == rank) {
1530: /* diagonal part A = c->A */
1531: ncols = ai[row-rstart+1] - ai[row-rstart];
1532: cols = aj + ai[row-rstart];
1533: if (!allcolumns) {
1534: for (k=0; k<ncols; k++) {
1535: #if defined(PETSC_USE_CTABLE)
1536: tcol = cmap_loc[cols[k]];
1537: #else
1538: tcol = cmap[cols[k]+cstart];
1539: #endif
1540: if (tcol) lens[j]++;
1541: }
1542: } else { /* allcolumns */
1543: lens[j] = ncols;
1544: }
1546: /* off-diagonal part B = c->B */
1547: ncols = bi[row-rstart+1] - bi[row-rstart];
1548: cols = bj + bi[row-rstart];
1549: if (!allcolumns) {
1550: for (k=0; k<ncols; k++) {
1551: #if defined(PETSC_USE_CTABLE)
1552: PetscTableFind(cmap,bmap[cols[k]]+1,&tcol);
1553: #else
1554: tcol = cmap[bmap[cols[k]]];
1555: #endif
1556: if (tcol) lens[j]++;
1557: }
1558: } else { /* allcolumns */
1559: lens[j] += ncols;
1560: }
1561: }
1562: }
1564: /* Create row map (rmap): global row of C -> local row of submat */
1565: #if defined(PETSC_USE_CTABLE)
1566: PetscTableCreate(nrow+1,C->rmap->N+1,&rmap);
1567: for (j=0; j<nrow; j++) {
1568: row = irow[j];
1569: proc = row2proc[j];
1570: if (proc == rank) { /* a local row */
1571: rmap_loc[row - rstart] = j;
1572: } else {
1573: PetscTableAdd(rmap,irow[j]+1,j+1,INSERT_VALUES);
1574: }
1575: }
1576: #else
1577: PetscCalloc1(C->rmap->N,&rmap);
1578: for (j=0; j<nrow; j++) {
1579: rmap[irow[j]] = j;
1580: }
1581: #endif
1583: /* Update lens from offproc data */
1584: /* recv a->j is done */
1585: MPI_Waitall(nrqs,r_waits3,r_status3);
1586: for (i=0; i<nrqs; i++) {
1587: proc = pa[i];
1588: sbuf1_i = sbuf1[proc];
1589: /* jmax = sbuf1_i[0]; if (jmax != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"jmax !=1"); */
1590: ct1 = 2 + 1;
1591: ct2 = 0;
1592: rbuf2_i = rbuf2[i]; /* received length of C->j */
1593: rbuf3_i = rbuf3[i]; /* received C->j */
1595: /* is_no = sbuf1_i[2*j-1]; if (is_no != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_no !=0"); */
1596: max1 = sbuf1_i[2];
1597: for (k=0; k<max1; k++,ct1++) {
1598: #if defined(PETSC_USE_CTABLE)
1599: PetscTableFind(rmap,sbuf1_i[ct1]+1,&row);
1600: row--;
1601: if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1602: #else
1603: row = rmap[sbuf1_i[ct1]]; /* the row index in submat */
1604: #endif
1605: /* Now, store row index of submat in sbuf1_i[ct1] */
1606: sbuf1_i[ct1] = row;
1608: nnz = rbuf2_i[ct1];
1609: if (!allcolumns) {
1610: for (l=0; l<nnz; l++,ct2++) {
1611: #if defined(PETSC_USE_CTABLE)
1612: if (rbuf3_i[ct2] >= cstart && rbuf3_i[ct2] <cend) {
1613: tcol = cmap_loc[rbuf3_i[ct2] - cstart];
1614: } else {
1615: PetscTableFind(cmap,rbuf3_i[ct2]+1,&tcol);
1616: }
1617: #else
1618: tcol = cmap[rbuf3_i[ct2]]; /* column index in submat */
1619: #endif
1620: if (tcol) lens[row]++;
1621: }
1622: } else { /* allcolumns */
1623: lens[row] += nnz;
1624: }
1625: }
1626: }
1627: MPI_Waitall(nrqr,s_waits3,s_status3);
1628: PetscFree4(r_waits3,s_waits3,r_status3,s_status3);
1630: /* Create the submatrices */
1631: MatCreate(PETSC_COMM_SELF,&submat);
1632: MatSetSizes(submat,nrow,ncol,PETSC_DETERMINE,PETSC_DETERMINE);
1634: ISGetBlockSize(isrow[0],&i);
1635: ISGetBlockSize(iscol[0],&j);
1636: MatSetBlockSizes(submat,i,j);
1637: MatSetType(submat,((PetscObject)A)->type_name);
1638: MatSeqAIJSetPreallocation(submat,0,lens);
1640: /* create struct Mat_SubSppt and attached it to submat */
1641: PetscNew(&smatis1);
1642: subc = (Mat_SeqAIJ*)submat->data;
1643: subc->submatis1 = smatis1;
1645: smatis1->id = 0;
1646: smatis1->nrqs = nrqs;
1647: smatis1->nrqr = nrqr;
1648: smatis1->rbuf1 = rbuf1;
1649: smatis1->rbuf2 = rbuf2;
1650: smatis1->rbuf3 = rbuf3;
1651: smatis1->sbuf2 = sbuf2;
1652: smatis1->req_source2 = req_source2;
1654: smatis1->sbuf1 = sbuf1;
1655: smatis1->ptr = ptr;
1656: smatis1->tmp = tmp;
1657: smatis1->ctr = ctr;
1659: smatis1->pa = pa;
1660: smatis1->req_size = req_size;
1661: smatis1->req_source1 = req_source1;
1663: smatis1->allcolumns = allcolumns;
1664: smatis1->singleis = PETSC_TRUE;
1665: smatis1->row2proc = row2proc;
1666: smatis1->rmap = rmap;
1667: smatis1->cmap = cmap;
1668: #if defined(PETSC_USE_CTABLE)
1669: smatis1->rmap_loc = rmap_loc;
1670: smatis1->cmap_loc = cmap_loc;
1671: #endif
1673: smatis1->destroy = submat->ops->destroy;
1674: submat->ops->destroy = MatDestroySubMatrix_SeqAIJ;
1675: submat->factortype = C->factortype;
1677: /* compute rmax */
1678: rmax = 0;
1679: for (i=0; i<nrow; i++) rmax = PetscMax(rmax,lens[i]);
1681: } else { /* scall == MAT_REUSE_MATRIX */
1682: submat = submats[0];
1683: if (submat->rmap->n != nrow || submat->cmap->n != ncol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1685: subc = (Mat_SeqAIJ*)submat->data;
1686: rmax = subc->rmax;
1687: smatis1 = subc->submatis1;
1688: nrqs = smatis1->nrqs;
1689: nrqr = smatis1->nrqr;
1690: rbuf1 = smatis1->rbuf1;
1691: rbuf2 = smatis1->rbuf2;
1692: rbuf3 = smatis1->rbuf3;
1693: req_source2 = smatis1->req_source2;
1695: sbuf1 = smatis1->sbuf1;
1696: sbuf2 = smatis1->sbuf2;
1697: ptr = smatis1->ptr;
1698: tmp = smatis1->tmp;
1699: ctr = smatis1->ctr;
1701: pa = smatis1->pa;
1702: req_size = smatis1->req_size;
1703: req_source1 = smatis1->req_source1;
1705: allcolumns = smatis1->allcolumns;
1706: row2proc = smatis1->row2proc;
1707: rmap = smatis1->rmap;
1708: cmap = smatis1->cmap;
1709: #if defined(PETSC_USE_CTABLE)
1710: rmap_loc = smatis1->rmap_loc;
1711: cmap_loc = smatis1->cmap_loc;
1712: #endif
1713: }
1715: /* Post recv matrix values */
1716: PetscMalloc3(nrqs+1,&rbuf4, rmax,&subcols, rmax,&subvals);
1717: PetscMalloc4(nrqs+1,&r_waits4,nrqr+1,&s_waits4,nrqs+1,&r_status4,nrqr+1,&s_status4);
1718: PetscObjectGetNewTag((PetscObject)C,&tag4);
1719: for (i=0; i<nrqs; ++i) {
1720: PetscMalloc1(rbuf2[i][0]+1,&rbuf4[i]);
1721: MPI_Irecv(rbuf4[i],rbuf2[i][0],MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);
1722: }
1724: /* Allocate sending buffers for a->a, and send them off */
1725: PetscMalloc1(nrqr+1,&sbuf_aa);
1726: for (i=0,j=0; i<nrqr; i++) j += req_size[i];
1727: PetscMalloc1(j+1,&sbuf_aa[0]);
1728: for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1];
1730: for (i=0; i<nrqr; i++) {
1731: rbuf1_i = rbuf1[i];
1732: sbuf_aa_i = sbuf_aa[i];
1733: ct1 = 2*rbuf1_i[0]+1;
1734: ct2 = 0;
1735: /* max1=rbuf1_i[0]; if (max1 != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"max1 !=1"); */
1737: kmax = rbuf1_i[2];
1738: for (k=0; k<kmax; k++,ct1++) {
1739: row = rbuf1_i[ct1] - rstart;
1740: nzA = ai[row+1] - ai[row];
1741: nzB = bi[row+1] - bi[row];
1742: ncols = nzA + nzB;
1743: cworkB = bj + bi[row];
1744: vworkA = a_a + ai[row];
1745: vworkB = b_a + bi[row];
1747: /* load the column values for this row into vals*/
1748: vals = sbuf_aa_i + ct2;
1750: lwrite = 0;
1751: for (l=0; l<nzB; l++) {
1752: if ((bmap[cworkB[l]]) < cstart) vals[lwrite++] = vworkB[l];
1753: }
1754: for (l=0; l<nzA; l++) vals[lwrite++] = vworkA[l];
1755: for (l=0; l<nzB; l++) {
1756: if ((bmap[cworkB[l]]) >= cend) vals[lwrite++] = vworkB[l];
1757: }
1759: ct2 += ncols;
1760: }
1761: MPI_Isend(sbuf_aa_i,req_size[i],MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);
1762: }
1764: /* Assemble submat */
1765: /* First assemble the local rows */
1766: for (j=0; j<nrow; j++) {
1767: row = irow[j];
1768: proc = row2proc[j];
1769: if (proc == rank) {
1770: Crow = row - rstart; /* local row index of C */
1771: #if defined(PETSC_USE_CTABLE)
1772: row = rmap_loc[Crow]; /* row index of submat */
1773: #else
1774: row = rmap[row];
1775: #endif
1777: if (allcolumns) {
1778: /* diagonal part A = c->A */
1779: ncols = ai[Crow+1] - ai[Crow];
1780: cols = aj + ai[Crow];
1781: vals = a->a + ai[Crow];
1782: i = 0;
1783: for (k=0; k<ncols; k++) {
1784: subcols[i] = cols[k] + cstart;
1785: subvals[i++] = vals[k];
1786: }
1788: /* off-diagonal part B = c->B */
1789: ncols = bi[Crow+1] - bi[Crow];
1790: cols = bj + bi[Crow];
1791: vals = b->a + bi[Crow];
1792: for (k=0; k<ncols; k++) {
1793: subcols[i] = bmap[cols[k]];
1794: subvals[i++] = vals[k];
1795: }
1797: MatSetValues_SeqAIJ(submat,1,&row,i,subcols,subvals,INSERT_VALUES);
1799: } else { /* !allcolumns */
1800: #if defined(PETSC_USE_CTABLE)
1801: /* diagonal part A = c->A */
1802: ncols = ai[Crow+1] - ai[Crow];
1803: cols = aj + ai[Crow];
1804: vals = a->a + ai[Crow];
1805: i = 0;
1806: for (k=0; k<ncols; k++) {
1807: tcol = cmap_loc[cols[k]];
1808: if (tcol) {
1809: subcols[i] = --tcol;
1810: subvals[i++] = vals[k];
1811: }
1812: }
1814: /* off-diagonal part B = c->B */
1815: ncols = bi[Crow+1] - bi[Crow];
1816: cols = bj + bi[Crow];
1817: vals = b->a + bi[Crow];
1818: for (k=0; k<ncols; k++) {
1819: PetscTableFind(cmap,bmap[cols[k]]+1,&tcol);
1820: if (tcol) {
1821: subcols[i] = --tcol;
1822: subvals[i++] = vals[k];
1823: }
1824: }
1825: #else
1826: /* diagonal part A = c->A */
1827: ncols = ai[Crow+1] - ai[Crow];
1828: cols = aj + ai[Crow];
1829: vals = a->a + ai[Crow];
1830: i = 0;
1831: for (k=0; k<ncols; k++) {
1832: tcol = cmap[cols[k]+cstart];
1833: if (tcol) {
1834: subcols[i] = --tcol;
1835: subvals[i++] = vals[k];
1836: }
1837: }
1839: /* off-diagonal part B = c->B */
1840: ncols = bi[Crow+1] - bi[Crow];
1841: cols = bj + bi[Crow];
1842: vals = b->a + bi[Crow];
1843: for (k=0; k<ncols; k++) {
1844: tcol = cmap[bmap[cols[k]]];
1845: if (tcol) {
1846: subcols[i] = --tcol;
1847: subvals[i++] = vals[k];
1848: }
1849: }
1850: #endif
1851: MatSetValues_SeqAIJ(submat,1,&row,i,subcols,subvals,INSERT_VALUES);
1852: }
1853: }
1854: }
1856: /* Now assemble the off-proc rows */
1857: for (i=0; i<nrqs; i++) { /* for each requested message */
1858: /* recv values from other processes */
1859: MPI_Waitany(nrqs,r_waits4,&idex,r_status4+i);
1860: proc = pa[idex];
1861: sbuf1_i = sbuf1[proc];
1862: /* jmax = sbuf1_i[0]; if (jmax != 1)SETERRQ1(PETSC_COMM_SELF,0,"jmax %d != 1",jmax); */
1863: ct1 = 2 + 1;
1864: ct2 = 0; /* count of received C->j */
1865: ct3 = 0; /* count of received C->j that will be inserted into submat */
1866: rbuf2_i = rbuf2[idex]; /* int** received length of C->j from other processes */
1867: rbuf3_i = rbuf3[idex]; /* int** received C->j from other processes */
1868: rbuf4_i = rbuf4[idex]; /* scalar** received C->a from other processes */
1870: /* is_no = sbuf1_i[2*j-1]; if (is_no != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_no !=0"); */
1871: max1 = sbuf1_i[2]; /* num of rows */
1872: for (k=0; k<max1; k++,ct1++) { /* for each recved row */
1873: row = sbuf1_i[ct1]; /* row index of submat */
1874: if (!allcolumns) {
1875: idex = 0;
1876: if (scall == MAT_INITIAL_MATRIX || !iscolsorted) {
1877: nnz = rbuf2_i[ct1]; /* num of C entries in this row */
1878: for (l=0; l<nnz; l++,ct2++) { /* for each recved column */
1879: #if defined(PETSC_USE_CTABLE)
1880: if (rbuf3_i[ct2] >= cstart && rbuf3_i[ct2] <cend) {
1881: tcol = cmap_loc[rbuf3_i[ct2] - cstart];
1882: } else {
1883: PetscTableFind(cmap,rbuf3_i[ct2]+1,&tcol);
1884: }
1885: #else
1886: tcol = cmap[rbuf3_i[ct2]];
1887: #endif
1888: if (tcol) {
1889: subcols[idex] = --tcol; /* may not be sorted */
1890: subvals[idex++] = rbuf4_i[ct2];
1892: /* We receive an entire column of C, but a subset of it needs to be inserted into submat.
1893: For reuse, we replace received C->j with index that should be inserted to submat */
1894: if (iscolsorted) rbuf3_i[ct3++] = ct2;
1895: }
1896: }
1897: MatSetValues_SeqAIJ(submat,1,&row,idex,subcols,subvals,INSERT_VALUES);
1898: } else { /* scall == MAT_REUSE_MATRIX */
1899: submat = submats[0];
1900: subc = (Mat_SeqAIJ*)submat->data;
1902: nnz = subc->i[row+1] - subc->i[row]; /* num of submat entries in this row */
1903: for (l=0; l<nnz; l++) {
1904: ct2 = rbuf3_i[ct3++]; /* index of rbuf4_i[] which needs to be inserted into submat */
1905: subvals[idex++] = rbuf4_i[ct2];
1906: }
1908: bj = subc->j + subc->i[row]; /* sorted column indices */
1909: MatSetValues_SeqAIJ(submat,1,&row,nnz,bj,subvals,INSERT_VALUES);
1910: }
1911: } else { /* allcolumns */
1912: nnz = rbuf2_i[ct1]; /* num of C entries in this row */
1913: MatSetValues_SeqAIJ(submat,1,&row,nnz,rbuf3_i+ct2,rbuf4_i+ct2,INSERT_VALUES);
1914: ct2 += nnz;
1915: }
1916: }
1917: }
1919: /* sending a->a are done */
1920: MPI_Waitall(nrqr,s_waits4,s_status4);
1921: PetscFree4(r_waits4,s_waits4,r_status4,s_status4);
1923: MatAssemblyBegin(submat,MAT_FINAL_ASSEMBLY);
1924: MatAssemblyEnd(submat,MAT_FINAL_ASSEMBLY);
1925: submats[0] = submat;
1927: /* Restore the indices */
1928: ISRestoreIndices(isrow[0],&irow);
1929: if (!allcolumns) {
1930: ISRestoreIndices(iscol[0],&icol);
1931: }
1933: /* Destroy allocated memory */
1934: for (i=0; i<nrqs; ++i) {
1935: PetscFree3(rbuf4[i],subcols,subvals);
1936: }
1937: PetscFree3(rbuf4,subcols,subvals);
1938: PetscFree(sbuf_aa[0]);
1939: PetscFree(sbuf_aa);
1941: if (scall == MAT_INITIAL_MATRIX) {
1942: PetscFree(lens);
1943: PetscFree(sbuf_aj[0]);
1944: PetscFree(sbuf_aj);
1945: }
1946: return(0);
1947: }
1949: PetscErrorCode MatCreateSubMatrices_MPIAIJ_SingleIS(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
1950: {
1952: PetscInt ncol;
1953: PetscBool colflag,allcolumns=PETSC_FALSE;
1956: /* Allocate memory to hold all the submatrices */
1957: if (scall == MAT_INITIAL_MATRIX) {
1958: PetscCalloc1(2,submat);
1959: }
1961: /* Check for special case: each processor gets entire matrix columns */
1962: ISIdentity(iscol[0],&colflag);
1963: ISGetLocalSize(iscol[0],&ncol);
1964: if (colflag && ncol == C->cmap->N) allcolumns = PETSC_TRUE;
1966: MatCreateSubMatrices_MPIAIJ_SingleIS_Local(C,ismax,isrow,iscol,scall,allcolumns,*submat);
1967: return(0);
1968: }
1970: PetscErrorCode MatCreateSubMatrices_MPIAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
1971: {
1973: PetscInt nmax,nstages=0,i,pos,max_no,nrow,ncol,in[2],out[2];
1974: PetscBool rowflag,colflag,wantallmatrix=PETSC_FALSE;
1975: Mat_SeqAIJ *subc;
1976: Mat_SubSppt *smat;
1979: /* Check for special case: each processor has a single IS */
1980: if (C->submat_singleis) { /* flag is set in PCSetUp_ASM() to skip MPIU_Allreduce() */
1981: MatCreateSubMatrices_MPIAIJ_SingleIS(C,ismax,isrow,iscol,scall,submat);
1982: C->submat_singleis = PETSC_FALSE; /* resume its default value in case C will be used for non-single IS */
1983: return(0);
1984: }
1986: /* Collect global wantallmatrix and nstages */
1987: if (!C->cmap->N) nmax=20*1000000/sizeof(PetscInt);
1988: else nmax = 20*1000000 / (C->cmap->N * sizeof(PetscInt));
1989: if (!nmax) nmax = 1;
1991: if (scall == MAT_INITIAL_MATRIX) {
1992: /* Collect global wantallmatrix and nstages */
1993: if (ismax == 1 && C->rmap->N == C->cmap->N) {
1994: ISIdentity(*isrow,&rowflag);
1995: ISIdentity(*iscol,&colflag);
1996: ISGetLocalSize(*isrow,&nrow);
1997: ISGetLocalSize(*iscol,&ncol);
1998: if (rowflag && colflag && nrow == C->rmap->N && ncol == C->cmap->N) {
1999: wantallmatrix = PETSC_TRUE;
2001: PetscOptionsGetBool(((PetscObject)C)->options,((PetscObject)C)->prefix,"-use_fast_submatrix",&wantallmatrix,NULL);
2002: }
2003: }
2005: /* Determine the number of stages through which submatrices are done
2006: Each stage will extract nmax submatrices.
2007: nmax is determined by the matrix column dimension.
2008: If the original matrix has 20M columns, only one submatrix per stage is allowed, etc.
2009: */
2010: nstages = ismax/nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */
2012: in[0] = -1*(PetscInt)wantallmatrix;
2013: in[1] = nstages;
2014: MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));
2015: wantallmatrix = (PetscBool)(-out[0]);
2016: nstages = out[1]; /* Make sure every processor loops through the global nstages */
2018: } else { /* MAT_REUSE_MATRIX */
2019: if (ismax) {
2020: subc = (Mat_SeqAIJ*)(*submat)[0]->data;
2021: smat = subc->submatis1;
2022: } else { /* (*submat)[0] is a dummy matrix */
2023: smat = (Mat_SubSppt*)(*submat)[0]->data;
2024: }
2025: if (!smat) {
2026: /* smat is not generated by MatCreateSubMatrix_MPIAIJ_All(...,MAT_INITIAL_MATRIX,...) */
2027: wantallmatrix = PETSC_TRUE;
2028: } else if (smat->singleis) {
2029: MatCreateSubMatrices_MPIAIJ_SingleIS(C,ismax,isrow,iscol,scall,submat);
2030: return(0);
2031: } else {
2032: nstages = smat->nstages;
2033: }
2034: }
2036: if (wantallmatrix) {
2037: MatCreateSubMatrix_MPIAIJ_All(C,MAT_GET_VALUES,scall,submat);
2038: return(0);
2039: }
2041: /* Allocate memory to hold all the submatrices and dummy submatrices */
2042: if (scall == MAT_INITIAL_MATRIX) {
2043: PetscCalloc1(ismax+nstages,submat);
2044: }
2046: for (i=0,pos=0; i<nstages; i++) {
2047: if (pos+nmax <= ismax) max_no = nmax;
2048: else if (pos >= ismax) max_no = 0;
2049: else max_no = ismax-pos;
2051: MatCreateSubMatrices_MPIAIJ_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);
2052: if (!max_no && scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */
2053: smat = (Mat_SubSppt*)(*submat)[pos]->data; pos++;
2054: smat->nstages = nstages;
2055: }
2056: pos += max_no;
2057: }
2059: if (ismax && scall == MAT_INITIAL_MATRIX) {
2060: /* save nstages for reuse */
2061: subc = (Mat_SeqAIJ*)(*submat)[0]->data;
2062: smat = subc->submatis1;
2063: smat->nstages = nstages;
2064: }
2065: return(0);
2066: }
2068: /* -------------------------------------------------------------------------*/
2069: PetscErrorCode MatCreateSubMatrices_MPIAIJ_Local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
2070: {
2071: Mat_MPIAIJ *c = (Mat_MPIAIJ*)C->data;
2072: Mat A = c->A;
2073: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)c->B->data,*subc;
2074: const PetscInt **icol,**irow;
2075: PetscInt *nrow,*ncol,start;
2077: PetscMPIInt rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr;
2078: PetscInt **sbuf1,**sbuf2,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1;
2079: PetscInt nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol;
2080: PetscInt **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2;
2081: PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
2082: #if defined(PETSC_USE_CTABLE)
2083: PetscTable *cmap,cmap_i=NULL,*rmap,rmap_i;
2084: #else
2085: PetscInt **cmap,*cmap_i=NULL,**rmap,*rmap_i;
2086: #endif
2087: const PetscInt *irow_i;
2088: PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
2089: MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
2090: MPI_Request *r_waits4,*s_waits3,*s_waits4;
2091: MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
2092: MPI_Status *r_status3,*r_status4,*s_status4;
2093: MPI_Comm comm;
2094: PetscScalar **rbuf4,*rbuf4_i,**sbuf_aa,*vals,*mat_a,*imat_a,*sbuf_aa_i;
2095: PetscMPIInt *onodes1,*olengths1,end;
2096: PetscInt **row2proc,*row2proc_i,ilen_row,*imat_ilen,*imat_j,*imat_i,old_row;
2097: Mat_SubSppt *smat_i;
2098: PetscBool *issorted,*allcolumns,colflag,iscsorted=PETSC_TRUE;
2099: PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen;
2102: PetscObjectGetComm((PetscObject)C,&comm);
2103: size = c->size;
2104: rank = c->rank;
2106: PetscMalloc4(ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax+1,&allcolumns);
2107: PetscMalloc5(ismax,(PetscInt***)&irow,ismax,(PetscInt***)&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);
2109: for (i=0; i<ismax; i++) {
2110: ISSorted(iscol[i],&issorted[i]);
2111: if (!issorted[i]) iscsorted = issorted[i];
2113: ISSorted(isrow[i],&issorted[i]);
2115: ISGetIndices(isrow[i],&irow[i]);
2116: ISGetLocalSize(isrow[i],&nrow[i]);
2118: /* Check for special case: allcolumn */
2119: ISIdentity(iscol[i],&colflag);
2120: ISGetLocalSize(iscol[i],&ncol[i]);
2121: if (colflag && ncol[i] == C->cmap->N) {
2122: allcolumns[i] = PETSC_TRUE;
2123: icol[i] = NULL;
2124: } else {
2125: allcolumns[i] = PETSC_FALSE;
2126: ISGetIndices(iscol[i],&icol[i]);
2127: }
2128: }
2130: if (scall == MAT_REUSE_MATRIX) {
2131: /* Assumes new rows are same length as the old rows */
2132: for (i=0; i<ismax; i++) {
2133: if (!submats[i]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"submats[%D] is null, cannot reuse",i);
2134: subc = (Mat_SeqAIJ*)submats[i]->data;
2135: 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");
2137: /* Initial matrix as if empty */
2138: PetscArrayzero(subc->ilen,submats[i]->rmap->n);
2140: smat_i = subc->submatis1;
2142: nrqs = smat_i->nrqs;
2143: nrqr = smat_i->nrqr;
2144: rbuf1 = smat_i->rbuf1;
2145: rbuf2 = smat_i->rbuf2;
2146: rbuf3 = smat_i->rbuf3;
2147: req_source2 = smat_i->req_source2;
2149: sbuf1 = smat_i->sbuf1;
2150: sbuf2 = smat_i->sbuf2;
2151: ptr = smat_i->ptr;
2152: tmp = smat_i->tmp;
2153: ctr = smat_i->ctr;
2155: pa = smat_i->pa;
2156: req_size = smat_i->req_size;
2157: req_source1 = smat_i->req_source1;
2159: allcolumns[i] = smat_i->allcolumns;
2160: row2proc[i] = smat_i->row2proc;
2161: rmap[i] = smat_i->rmap;
2162: cmap[i] = smat_i->cmap;
2163: }
2165: if (!ismax){ /* Get dummy submatrices and retrieve struct submatis1 */
2166: if (!submats[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"submats are null, cannot reuse");
2167: smat_i = (Mat_SubSppt*)submats[0]->data;
2169: nrqs = smat_i->nrqs;
2170: nrqr = smat_i->nrqr;
2171: rbuf1 = smat_i->rbuf1;
2172: rbuf2 = smat_i->rbuf2;
2173: rbuf3 = smat_i->rbuf3;
2174: req_source2 = smat_i->req_source2;
2176: sbuf1 = smat_i->sbuf1;
2177: sbuf2 = smat_i->sbuf2;
2178: ptr = smat_i->ptr;
2179: tmp = smat_i->tmp;
2180: ctr = smat_i->ctr;
2182: pa = smat_i->pa;
2183: req_size = smat_i->req_size;
2184: req_source1 = smat_i->req_source1;
2186: allcolumns[0] = PETSC_FALSE;
2187: }
2188: } else { /* scall == MAT_INITIAL_MATRIX */
2189: /* Get some new tags to keep the communication clean */
2190: PetscObjectGetNewTag((PetscObject)C,&tag2);
2191: PetscObjectGetNewTag((PetscObject)C,&tag3);
2193: /* evaluate communication - mesg to who, length of mesg, and buffer space
2194: required. Based on this, buffers are allocated, and data copied into them*/
2195: PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4); /* mesg size, initialize work vectors */
2197: for (i=0; i<ismax; i++) {
2198: jmax = nrow[i];
2199: irow_i = irow[i];
2201: PetscMalloc1(jmax,&row2proc_i);
2202: row2proc[i] = row2proc_i;
2204: if (issorted[i]) proc = 0;
2205: for (j=0; j<jmax; j++) {
2206: if (!issorted[i]) proc = 0;
2207: row = irow_i[j];
2208: while (row >= C->rmap->range[proc+1]) proc++;
2209: w4[proc]++;
2210: row2proc_i[j] = proc; /* map row index to proc */
2211: }
2212: for (j=0; j<size; j++) {
2213: if (w4[j]) { w1[j] += w4[j]; w3[j]++; w4[j] = 0;}
2214: }
2215: }
2217: nrqs = 0; /* no of outgoing messages */
2218: msz = 0; /* total mesg length (for all procs) */
2219: w1[rank] = 0; /* no mesg sent to self */
2220: w3[rank] = 0;
2221: for (i=0; i<size; i++) {
2222: if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
2223: }
2224: PetscMalloc1(nrqs+1,&pa); /*(proc -array)*/
2225: for (i=0,j=0; i<size; i++) {
2226: if (w1[i]) { pa[j] = i; j++; }
2227: }
2229: /* Each message would have a header = 1 + 2*(no of IS) + data */
2230: for (i=0; i<nrqs; i++) {
2231: j = pa[i];
2232: w1[j] += w2[j] + 2* w3[j];
2233: msz += w1[j];
2234: }
2235: PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);
2237: /* Determine the number of messages to expect, their lengths, from from-ids */
2238: PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);
2239: PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);
2241: /* Now post the Irecvs corresponding to these messages */
2242: tag0 = ((PetscObject)C)->tag;
2243: PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);
2245: PetscFree(onodes1);
2246: PetscFree(olengths1);
2248: /* Allocate Memory for outgoing messages */
2249: PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);
2250: PetscArrayzero(sbuf1,size);
2251: PetscArrayzero(ptr,size);
2253: {
2254: PetscInt *iptr = tmp;
2255: k = 0;
2256: for (i=0; i<nrqs; i++) {
2257: j = pa[i];
2258: iptr += k;
2259: sbuf1[j] = iptr;
2260: k = w1[j];
2261: }
2262: }
2264: /* Form the outgoing messages. Initialize the header space */
2265: for (i=0; i<nrqs; i++) {
2266: j = pa[i];
2267: sbuf1[j][0] = 0;
2268: PetscArrayzero(sbuf1[j]+1,2*w3[j]);
2269: ptr[j] = sbuf1[j] + 2*w3[j] + 1;
2270: }
2272: /* Parse the isrow and copy data into outbuf */
2273: for (i=0; i<ismax; i++) {
2274: row2proc_i = row2proc[i];
2275: PetscArrayzero(ctr,size);
2276: irow_i = irow[i];
2277: jmax = nrow[i];
2278: for (j=0; j<jmax; j++) { /* parse the indices of each IS */
2279: proc = row2proc_i[j];
2280: if (proc != rank) { /* copy to the outgoing buf*/
2281: ctr[proc]++;
2282: *ptr[proc] = irow_i[j];
2283: ptr[proc]++;
2284: }
2285: }
2286: /* Update the headers for the current IS */
2287: for (j=0; j<size; j++) { /* Can Optimise this loop too */
2288: if ((ctr_j = ctr[j])) {
2289: sbuf1_j = sbuf1[j];
2290: k = ++sbuf1_j[0];
2291: sbuf1_j[2*k] = ctr_j;
2292: sbuf1_j[2*k-1] = i;
2293: }
2294: }
2295: }
2297: /* Now post the sends */
2298: PetscMalloc1(nrqs+1,&s_waits1);
2299: for (i=0; i<nrqs; ++i) {
2300: j = pa[i];
2301: MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);
2302: }
2304: /* Post Receives to capture the buffer size */
2305: PetscMalloc1(nrqs+1,&r_waits2);
2306: PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);
2307: rbuf2[0] = tmp + msz;
2308: for (i=1; i<nrqs; ++i) {
2309: rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
2310: }
2311: for (i=0; i<nrqs; ++i) {
2312: j = pa[i];
2313: MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);
2314: }
2316: /* Send to other procs the buf size they should allocate */
2317: /* Receive messages*/
2318: PetscMalloc1(nrqr+1,&s_waits2);
2319: PetscMalloc1(nrqr+1,&r_status1);
2320: PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);
2321: {
2322: PetscInt *sAi = a->i,*sBi = b->i,id,rstart = C->rmap->rstart;
2323: PetscInt *sbuf2_i;
2325: MPI_Waitall(nrqr,r_waits1,r_status1);
2326: for (i=0; i<nrqr; ++i) {
2327: req_size[i] = 0;
2328: rbuf1_i = rbuf1[i];
2329: start = 2*rbuf1_i[0] + 1;
2330: MPI_Get_count(r_status1+i,MPIU_INT,&end);
2331: PetscMalloc1(end+1,&sbuf2[i]);
2332: sbuf2_i = sbuf2[i];
2333: for (j=start; j<end; j++) {
2334: id = rbuf1_i[j] - rstart;
2335: ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
2336: sbuf2_i[j] = ncols;
2337: req_size[i] += ncols;
2338: }
2339: req_source1[i] = r_status1[i].MPI_SOURCE;
2340: /* form the header */
2341: sbuf2_i[0] = req_size[i];
2342: for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];
2344: MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);
2345: }
2346: }
2347: PetscFree(r_status1);
2348: PetscFree(r_waits1);
2349: PetscFree4(w1,w2,w3,w4);
2351: /* Receive messages*/
2352: PetscMalloc1(nrqs+1,&r_waits3);
2353: PetscMalloc1(nrqs+1,&r_status2);
2355: MPI_Waitall(nrqs,r_waits2,r_status2);
2356: for (i=0; i<nrqs; ++i) {
2357: PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);
2358: req_source2[i] = r_status2[i].MPI_SOURCE;
2359: MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);
2360: }
2361: PetscFree(r_status2);
2362: PetscFree(r_waits2);
2364: /* Wait on sends1 and sends2 */
2365: PetscMalloc1(nrqs+1,&s_status1);
2366: PetscMalloc1(nrqr+1,&s_status2);
2368: if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status1);}
2369: if (nrqr) {MPI_Waitall(nrqr,s_waits2,s_status2);}
2370: PetscFree(s_status1);
2371: PetscFree(s_status2);
2372: PetscFree(s_waits1);
2373: PetscFree(s_waits2);
2375: /* Now allocate sending buffers for a->j, and send them off */
2376: PetscMalloc1(nrqr+1,&sbuf_aj);
2377: for (i=0,j=0; i<nrqr; i++) j += req_size[i];
2378: PetscMalloc1(j+1,&sbuf_aj[0]);
2379: for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
2381: PetscMalloc1(nrqr+1,&s_waits3);
2382: {
2383: PetscInt nzA,nzB,*a_i = a->i,*b_i = b->i,lwrite;
2384: PetscInt *cworkA,*cworkB,cstart = C->cmap->rstart,rstart = C->rmap->rstart,*bmap = c->garray;
2385: PetscInt cend = C->cmap->rend;
2386: PetscInt *a_j = a->j,*b_j = b->j,ctmp;
2388: for (i=0; i<nrqr; i++) {
2389: rbuf1_i = rbuf1[i];
2390: sbuf_aj_i = sbuf_aj[i];
2391: ct1 = 2*rbuf1_i[0] + 1;
2392: ct2 = 0;
2393: for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
2394: kmax = rbuf1[i][2*j];
2395: for (k=0; k<kmax; k++,ct1++) {
2396: row = rbuf1_i[ct1] - rstart;
2397: nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row];
2398: ncols = nzA + nzB;
2399: cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];
2401: /* load the column indices for this row into cols */
2402: cols = sbuf_aj_i + ct2;
2404: lwrite = 0;
2405: for (l=0; l<nzB; l++) {
2406: if ((ctmp = bmap[cworkB[l]]) < cstart) cols[lwrite++] = ctmp;
2407: }
2408: for (l=0; l<nzA; l++) cols[lwrite++] = cstart + cworkA[l];
2409: for (l=0; l<nzB; l++) {
2410: if ((ctmp = bmap[cworkB[l]]) >= cend) cols[lwrite++] = ctmp;
2411: }
2413: ct2 += ncols;
2414: }
2415: }
2416: MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);
2417: }
2418: }
2419: PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);
2421: /* create col map: global col of C -> local col of submatrices */
2422: {
2423: const PetscInt *icol_i;
2424: #if defined(PETSC_USE_CTABLE)
2425: for (i=0; i<ismax; i++) {
2426: if (!allcolumns[i]) {
2427: PetscTableCreate(ncol[i]+1,C->cmap->N+1,&cmap[i]);
2429: jmax = ncol[i];
2430: icol_i = icol[i];
2431: cmap_i = cmap[i];
2432: for (j=0; j<jmax; j++) {
2433: PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);
2434: }
2435: } else cmap[i] = NULL;
2436: }
2437: #else
2438: for (i=0; i<ismax; i++) {
2439: if (!allcolumns[i]) {
2440: PetscCalloc1(C->cmap->N,&cmap[i]);
2441: jmax = ncol[i];
2442: icol_i = icol[i];
2443: cmap_i = cmap[i];
2444: for (j=0; j<jmax; j++) {
2445: cmap_i[icol_i[j]] = j+1;
2446: }
2447: } else cmap[i] = NULL;
2448: }
2449: #endif
2450: }
2452: /* Create lens which is required for MatCreate... */
2453: for (i=0,j=0; i<ismax; i++) j += nrow[i];
2454: PetscMalloc1(ismax,&lens);
2456: if (ismax) {
2457: PetscCalloc1(j,&lens[0]);
2458: }
2459: for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1];
2461: /* Update lens from local data */
2462: for (i=0; i<ismax; i++) {
2463: row2proc_i = row2proc[i];
2464: jmax = nrow[i];
2465: if (!allcolumns[i]) cmap_i = cmap[i];
2466: irow_i = irow[i];
2467: lens_i = lens[i];
2468: for (j=0; j<jmax; j++) {
2469: row = irow_i[j];
2470: proc = row2proc_i[j];
2471: if (proc == rank) {
2472: MatGetRow_MPIAIJ(C,row,&ncols,&cols,0);
2473: if (!allcolumns[i]) {
2474: for (k=0; k<ncols; k++) {
2475: #if defined(PETSC_USE_CTABLE)
2476: PetscTableFind(cmap_i,cols[k]+1,&tcol);
2477: #else
2478: tcol = cmap_i[cols[k]];
2479: #endif
2480: if (tcol) lens_i[j]++;
2481: }
2482: } else { /* allcolumns */
2483: lens_i[j] = ncols;
2484: }
2485: MatRestoreRow_MPIAIJ(C,row,&ncols,&cols,0);
2486: }
2487: }
2488: }
2490: /* Create row map: global row of C -> local row of submatrices */
2491: #if defined(PETSC_USE_CTABLE)
2492: for (i=0; i<ismax; i++) {
2493: PetscTableCreate(nrow[i]+1,C->rmap->N+1,&rmap[i]);
2494: irow_i = irow[i];
2495: jmax = nrow[i];
2496: for (j=0; j<jmax; j++) {
2497: PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);
2498: }
2499: }
2500: #else
2501: for (i=0; i<ismax; i++) {
2502: PetscCalloc1(C->rmap->N,&rmap[i]);
2503: rmap_i = rmap[i];
2504: irow_i = irow[i];
2505: jmax = nrow[i];
2506: for (j=0; j<jmax; j++) {
2507: rmap_i[irow_i[j]] = j;
2508: }
2509: }
2510: #endif
2512: /* Update lens from offproc data */
2513: {
2514: PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i;
2516: MPI_Waitall(nrqs,r_waits3,r_status3);
2517: for (tmp2=0; tmp2<nrqs; tmp2++) {
2518: sbuf1_i = sbuf1[pa[tmp2]];
2519: jmax = sbuf1_i[0];
2520: ct1 = 2*jmax+1;
2521: ct2 = 0;
2522: rbuf2_i = rbuf2[tmp2];
2523: rbuf3_i = rbuf3[tmp2];
2524: for (j=1; j<=jmax; j++) {
2525: is_no = sbuf1_i[2*j-1];
2526: max1 = sbuf1_i[2*j];
2527: lens_i = lens[is_no];
2528: if (!allcolumns[is_no]) cmap_i = cmap[is_no];
2529: rmap_i = rmap[is_no];
2530: for (k=0; k<max1; k++,ct1++) {
2531: #if defined(PETSC_USE_CTABLE)
2532: PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);
2533: row--;
2534: if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
2535: #else
2536: row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
2537: #endif
2538: max2 = rbuf2_i[ct1];
2539: for (l=0; l<max2; l++,ct2++) {
2540: if (!allcolumns[is_no]) {
2541: #if defined(PETSC_USE_CTABLE)
2542: PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);
2543: #else
2544: tcol = cmap_i[rbuf3_i[ct2]];
2545: #endif
2546: if (tcol) lens_i[row]++;
2547: } else { /* allcolumns */
2548: lens_i[row]++; /* lens_i[row] += max2 ? */
2549: }
2550: }
2551: }
2552: }
2553: }
2554: }
2555: PetscFree(r_waits3);
2556: if (nrqr) {MPI_Waitall(nrqr,s_waits3,s_status3);}
2557: PetscFree2(r_status3,s_status3);
2558: PetscFree(s_waits3);
2560: /* Create the submatrices */
2561: for (i=0; i<ismax; i++) {
2562: PetscInt rbs,cbs;
2564: ISGetBlockSize(isrow[i],&rbs);
2565: ISGetBlockSize(iscol[i],&cbs);
2567: MatCreate(PETSC_COMM_SELF,submats+i);
2568: MatSetSizes(submats[i],nrow[i],ncol[i],PETSC_DETERMINE,PETSC_DETERMINE);
2570: MatSetBlockSizes(submats[i],rbs,cbs);
2571: MatSetType(submats[i],((PetscObject)A)->type_name);
2572: MatSeqAIJSetPreallocation(submats[i],0,lens[i]);
2574: /* create struct Mat_SubSppt and attached it to submat */
2575: PetscNew(&smat_i);
2576: subc = (Mat_SeqAIJ*)submats[i]->data;
2577: subc->submatis1 = smat_i;
2579: smat_i->destroy = submats[i]->ops->destroy;
2580: submats[i]->ops->destroy = MatDestroySubMatrix_SeqAIJ;
2581: submats[i]->factortype = C->factortype;
2583: smat_i->id = i;
2584: smat_i->nrqs = nrqs;
2585: smat_i->nrqr = nrqr;
2586: smat_i->rbuf1 = rbuf1;
2587: smat_i->rbuf2 = rbuf2;
2588: smat_i->rbuf3 = rbuf3;
2589: smat_i->sbuf2 = sbuf2;
2590: smat_i->req_source2 = req_source2;
2592: smat_i->sbuf1 = sbuf1;
2593: smat_i->ptr = ptr;
2594: smat_i->tmp = tmp;
2595: smat_i->ctr = ctr;
2597: smat_i->pa = pa;
2598: smat_i->req_size = req_size;
2599: smat_i->req_source1 = req_source1;
2601: smat_i->allcolumns = allcolumns[i];
2602: smat_i->singleis = PETSC_FALSE;
2603: smat_i->row2proc = row2proc[i];
2604: smat_i->rmap = rmap[i];
2605: smat_i->cmap = cmap[i];
2606: }
2608: if (!ismax) { /* Create dummy submats[0] for reuse struct subc */
2609: MatCreate(PETSC_COMM_SELF,&submats[0]);
2610: MatSetSizes(submats[0],0,0,PETSC_DETERMINE,PETSC_DETERMINE);
2611: MatSetType(submats[0],MATDUMMY);
2613: /* create struct Mat_SubSppt and attached it to submat */
2614: PetscNewLog(submats[0],&smat_i);
2615: submats[0]->data = (void*)smat_i;
2617: smat_i->destroy = submats[0]->ops->destroy;
2618: submats[0]->ops->destroy = MatDestroySubMatrix_Dummy;
2619: submats[0]->factortype = C->factortype;
2621: smat_i->id = 0;
2622: smat_i->nrqs = nrqs;
2623: smat_i->nrqr = nrqr;
2624: smat_i->rbuf1 = rbuf1;
2625: smat_i->rbuf2 = rbuf2;
2626: smat_i->rbuf3 = rbuf3;
2627: smat_i->sbuf2 = sbuf2;
2628: smat_i->req_source2 = req_source2;
2630: smat_i->sbuf1 = sbuf1;
2631: smat_i->ptr = ptr;
2632: smat_i->tmp = tmp;
2633: smat_i->ctr = ctr;
2635: smat_i->pa = pa;
2636: smat_i->req_size = req_size;
2637: smat_i->req_source1 = req_source1;
2639: smat_i->allcolumns = PETSC_FALSE;
2640: smat_i->singleis = PETSC_FALSE;
2641: smat_i->row2proc = NULL;
2642: smat_i->rmap = NULL;
2643: smat_i->cmap = NULL;
2644: }
2646: if (ismax) {PetscFree(lens[0]);}
2647: PetscFree(lens);
2648: PetscFree(sbuf_aj[0]);
2649: PetscFree(sbuf_aj);
2651: } /* endof scall == MAT_INITIAL_MATRIX */
2653: /* Post recv matrix values */
2654: PetscObjectGetNewTag((PetscObject)C,&tag4);
2655: PetscMalloc1(nrqs+1,&rbuf4);
2656: PetscMalloc1(nrqs+1,&r_waits4);
2657: PetscMalloc1(nrqs+1,&r_status4);
2658: PetscMalloc1(nrqr+1,&s_status4);
2659: for (i=0; i<nrqs; ++i) {
2660: PetscMalloc1(rbuf2[i][0]+1,&rbuf4[i]);
2661: MPI_Irecv(rbuf4[i],rbuf2[i][0],MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);
2662: }
2664: /* Allocate sending buffers for a->a, and send them off */
2665: PetscMalloc1(nrqr+1,&sbuf_aa);
2666: for (i=0,j=0; i<nrqr; i++) j += req_size[i];
2667: PetscMalloc1(j+1,&sbuf_aa[0]);
2668: for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1];
2670: PetscMalloc1(nrqr+1,&s_waits4);
2671: {
2672: PetscInt nzA,nzB,*a_i = a->i,*b_i = b->i, *cworkB,lwrite;
2673: PetscInt cstart = C->cmap->rstart,rstart = C->rmap->rstart,*bmap = c->garray;
2674: PetscInt cend = C->cmap->rend;
2675: PetscInt *b_j = b->j;
2676: PetscScalar *vworkA,*vworkB,*a_a = a->a,*b_a = b->a;
2678: for (i=0; i<nrqr; i++) {
2679: rbuf1_i = rbuf1[i];
2680: sbuf_aa_i = sbuf_aa[i];
2681: ct1 = 2*rbuf1_i[0]+1;
2682: ct2 = 0;
2683: for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
2684: kmax = rbuf1_i[2*j];
2685: for (k=0; k<kmax; k++,ct1++) {
2686: row = rbuf1_i[ct1] - rstart;
2687: nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row];
2688: ncols = nzA + nzB;
2689: cworkB = b_j + b_i[row];
2690: vworkA = a_a + a_i[row];
2691: vworkB = b_a + b_i[row];
2693: /* load the column values for this row into vals*/
2694: vals = sbuf_aa_i+ct2;
2696: lwrite = 0;
2697: for (l=0; l<nzB; l++) {
2698: if ((bmap[cworkB[l]]) < cstart) vals[lwrite++] = vworkB[l];
2699: }
2700: for (l=0; l<nzA; l++) vals[lwrite++] = vworkA[l];
2701: for (l=0; l<nzB; l++) {
2702: if ((bmap[cworkB[l]]) >= cend) vals[lwrite++] = vworkB[l];
2703: }
2705: ct2 += ncols;
2706: }
2707: }
2708: MPI_Isend(sbuf_aa_i,req_size[i],MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);
2709: }
2710: }
2712: /* Assemble the matrices */
2713: /* First assemble the local rows */
2714: for (i=0; i<ismax; i++) {
2715: row2proc_i = row2proc[i];
2716: subc = (Mat_SeqAIJ*)submats[i]->data;
2717: imat_ilen = subc->ilen;
2718: imat_j = subc->j;
2719: imat_i = subc->i;
2720: imat_a = subc->a;
2722: if (!allcolumns[i]) cmap_i = cmap[i];
2723: rmap_i = rmap[i];
2724: irow_i = irow[i];
2725: jmax = nrow[i];
2726: for (j=0; j<jmax; j++) {
2727: row = irow_i[j];
2728: proc = row2proc_i[j];
2729: if (proc == rank) {
2730: old_row = row;
2731: #if defined(PETSC_USE_CTABLE)
2732: PetscTableFind(rmap_i,row+1,&row);
2733: row--;
2734: #else
2735: row = rmap_i[row];
2736: #endif
2737: ilen_row = imat_ilen[row];
2738: MatGetRow_MPIAIJ(C,old_row,&ncols,&cols,&vals);
2739: mat_i = imat_i[row];
2740: mat_a = imat_a + mat_i;
2741: mat_j = imat_j + mat_i;
2742: if (!allcolumns[i]) {
2743: for (k=0; k<ncols; k++) {
2744: #if defined(PETSC_USE_CTABLE)
2745: PetscTableFind(cmap_i,cols[k]+1,&tcol);
2746: #else
2747: tcol = cmap_i[cols[k]];
2748: #endif
2749: if (tcol) {
2750: *mat_j++ = tcol - 1;
2751: *mat_a++ = vals[k];
2752: ilen_row++;
2753: }
2754: }
2755: } else { /* allcolumns */
2756: for (k=0; k<ncols; k++) {
2757: *mat_j++ = cols[k]; /* global col index! */
2758: *mat_a++ = vals[k];
2759: ilen_row++;
2760: }
2761: }
2762: MatRestoreRow_MPIAIJ(C,old_row,&ncols,&cols,&vals);
2764: imat_ilen[row] = ilen_row;
2765: }
2766: }
2767: }
2769: /* Now assemble the off proc rows */
2770: MPI_Waitall(nrqs,r_waits4,r_status4);
2771: for (tmp2=0; tmp2<nrqs; tmp2++) {
2772: sbuf1_i = sbuf1[pa[tmp2]];
2773: jmax = sbuf1_i[0];
2774: ct1 = 2*jmax + 1;
2775: ct2 = 0;
2776: rbuf2_i = rbuf2[tmp2];
2777: rbuf3_i = rbuf3[tmp2];
2778: rbuf4_i = rbuf4[tmp2];
2779: for (j=1; j<=jmax; j++) {
2780: is_no = sbuf1_i[2*j-1];
2781: rmap_i = rmap[is_no];
2782: if (!allcolumns[is_no]) cmap_i = cmap[is_no];
2783: subc = (Mat_SeqAIJ*)submats[is_no]->data;
2784: imat_ilen = subc->ilen;
2785: imat_j = subc->j;
2786: imat_i = subc->i;
2787: imat_a = subc->a;
2788: max1 = sbuf1_i[2*j];
2789: for (k=0; k<max1; k++,ct1++) {
2790: row = sbuf1_i[ct1];
2791: #if defined(PETSC_USE_CTABLE)
2792: PetscTableFind(rmap_i,row+1,&row);
2793: row--;
2794: #else
2795: row = rmap_i[row];
2796: #endif
2797: ilen = imat_ilen[row];
2798: mat_i = imat_i[row];
2799: mat_a = imat_a + mat_i;
2800: mat_j = imat_j + mat_i;
2801: max2 = rbuf2_i[ct1];
2802: if (!allcolumns[is_no]) {
2803: for (l=0; l<max2; l++,ct2++) {
2804: #if defined(PETSC_USE_CTABLE)
2805: PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);
2806: #else
2807: tcol = cmap_i[rbuf3_i[ct2]];
2808: #endif
2809: if (tcol) {
2810: *mat_j++ = tcol - 1;
2811: *mat_a++ = rbuf4_i[ct2];
2812: ilen++;
2813: }
2814: }
2815: } else { /* allcolumns */
2816: for (l=0; l<max2; l++,ct2++) {
2817: *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */
2818: *mat_a++ = rbuf4_i[ct2];
2819: ilen++;
2820: }
2821: }
2822: imat_ilen[row] = ilen;
2823: }
2824: }
2825: }
2827: if (!iscsorted) { /* sort column indices of the rows */
2828: for (i=0; i<ismax; i++) {
2829: subc = (Mat_SeqAIJ*)submats[i]->data;
2830: imat_j = subc->j;
2831: imat_i = subc->i;
2832: imat_a = subc->a;
2833: imat_ilen = subc->ilen;
2835: if (allcolumns[i]) continue;
2836: jmax = nrow[i];
2837: for (j=0; j<jmax; j++) {
2838: mat_i = imat_i[j];
2839: mat_a = imat_a + mat_i;
2840: mat_j = imat_j + mat_i;
2841: PetscSortIntWithScalarArray(imat_ilen[j],mat_j,mat_a);
2842: }
2843: }
2844: }
2846: PetscFree(r_status4);
2847: PetscFree(r_waits4);
2848: if (nrqr) {MPI_Waitall(nrqr,s_waits4,s_status4);}
2849: PetscFree(s_waits4);
2850: PetscFree(s_status4);
2852: /* Restore the indices */
2853: for (i=0; i<ismax; i++) {
2854: ISRestoreIndices(isrow[i],irow+i);
2855: if (!allcolumns[i]) {
2856: ISRestoreIndices(iscol[i],icol+i);
2857: }
2858: }
2860: for (i=0; i<ismax; i++) {
2861: MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);
2862: MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);
2863: }
2865: /* Destroy allocated memory */
2866: PetscFree(sbuf_aa[0]);
2867: PetscFree(sbuf_aa);
2868: PetscFree5(*(PetscInt***)&irow,*(PetscInt***)&icol,nrow,ncol,issorted);
2870: for (i=0; i<nrqs; ++i) {
2871: PetscFree(rbuf4[i]);
2872: }
2873: PetscFree(rbuf4);
2875: PetscFree4(row2proc,cmap,rmap,allcolumns);
2876: return(0);
2877: }
2879: /*
2880: Permute A & B into C's *local* index space using rowemb,dcolemb for A and rowemb,ocolemb for B.
2881: Embeddings are supposed to be injections and the above implies that the range of rowemb is a subset
2882: of [0,m), dcolemb is in [0,n) and ocolemb is in [N-n).
2883: If pattern == DIFFERENT_NONZERO_PATTERN, C is preallocated according to A&B.
2884: After that B's columns are mapped into C's global column space, so that C is in the "disassembled"
2885: state, and needs to be "assembled" later by compressing B's column space.
2887: This function may be called in lieu of preallocation, so C should not be expected to be preallocated.
2888: Following this call, C->A & C->B have been created, even if empty.
2889: */
2890: PetscErrorCode MatSetSeqMats_MPIAIJ(Mat C,IS rowemb,IS dcolemb,IS ocolemb,MatStructure pattern,Mat A,Mat B)
2891: {
2892: /* If making this function public, change the error returned in this function away from _PLIB. */
2894: Mat_MPIAIJ *aij;
2895: Mat_SeqAIJ *Baij;
2896: PetscBool seqaij,Bdisassembled;
2897: PetscInt m,n,*nz,i,j,ngcol,col,rstart,rend,shift,count;
2898: PetscScalar v;
2899: const PetscInt *rowindices,*colindices;
2902: /* Check to make sure the component matrices (and embeddings) are compatible with C. */
2903: if (A) {
2904: PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&seqaij);
2905: if (!seqaij) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Diagonal matrix is of wrong type");
2906: if (rowemb) {
2907: ISGetLocalSize(rowemb,&m);
2908: 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);
2909: } else {
2910: if (C->rmap->n != A->rmap->n) {
2911: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Diag seq matrix is row-incompatible with the MPIAIJ matrix");
2912: }
2913: }
2914: if (dcolemb) {
2915: ISGetLocalSize(dcolemb,&n);
2916: 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);
2917: } else {
2918: if (C->cmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Diag seq matrix is col-incompatible with the MPIAIJ matrix");
2919: }
2920: }
2921: if (B) {
2922: PetscObjectBaseTypeCompare((PetscObject)B,MATSEQAIJ,&seqaij);
2923: if (!seqaij) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Off-diagonal matrix is of wrong type");
2924: if (rowemb) {
2925: ISGetLocalSize(rowemb,&m);
2926: 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);
2927: } else {
2928: if (C->rmap->n != B->rmap->n) {
2929: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Off-diag seq matrix is row-incompatible with the MPIAIJ matrix");
2930: }
2931: }
2932: if (ocolemb) {
2933: ISGetLocalSize(ocolemb,&n);
2934: 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);
2935: } else {
2936: 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");
2937: }
2938: }
2940: aij = (Mat_MPIAIJ*)C->data;
2941: if (!aij->A) {
2942: /* Mimic parts of MatMPIAIJSetPreallocation() */
2943: MatCreate(PETSC_COMM_SELF,&aij->A);
2944: MatSetSizes(aij->A,C->rmap->n,C->cmap->n,C->rmap->n,C->cmap->n);
2945: MatSetBlockSizesFromMats(aij->A,C,C);
2946: MatSetType(aij->A,MATSEQAIJ);
2947: PetscLogObjectParent((PetscObject)C,(PetscObject)aij->A);
2948: }
2949: if (A) {
2950: MatSetSeqMat_SeqAIJ(aij->A,rowemb,dcolemb,pattern,A);
2951: } else {
2952: MatSetUp(aij->A);
2953: }
2954: if (B) { /* Destroy the old matrix or the column map, depending on the sparsity pattern. */
2955: /*
2956: If pattern == DIFFERENT_NONZERO_PATTERN, we reallocate B and
2957: need to "disassemble" B -- convert it to using C's global indices.
2958: To insert the values we take the safer, albeit more expensive, route of MatSetValues().
2960: If pattern == SUBSET_NONZERO_PATTERN, we do not "disassemble" B and do not reallocate;
2961: we MatZeroValues(B) first, so there may be a bunch of zeros that, perhaps, could be compacted out.
2963: TODO: Put B's values into aij->B's aij structure in place using the embedding ISs?
2964: At least avoid calling MatSetValues() and the implied searches?
2965: */
2967: if (B && pattern == DIFFERENT_NONZERO_PATTERN) {
2968: #if defined(PETSC_USE_CTABLE)
2969: PetscTableDestroy(&aij->colmap);
2970: #else
2971: PetscFree(aij->colmap);
2972: /* A bit of a HACK: ideally we should deal with case aij->B all in one code block below. */
2973: if (aij->B) {
2974: PetscLogObjectMemory((PetscObject)C,-aij->B->cmap->n*sizeof(PetscInt));
2975: }
2976: #endif
2977: ngcol = 0;
2978: if (aij->lvec) {
2979: VecGetSize(aij->lvec,&ngcol);
2980: }
2981: if (aij->garray) {
2982: PetscFree(aij->garray);
2983: PetscLogObjectMemory((PetscObject)C,-ngcol*sizeof(PetscInt));
2984: }
2985: VecDestroy(&aij->lvec);
2986: VecScatterDestroy(&aij->Mvctx);
2987: }
2988: if (aij->B && B && pattern == DIFFERENT_NONZERO_PATTERN) {
2989: MatDestroy(&aij->B);
2990: }
2991: if (aij->B && B && pattern == SUBSET_NONZERO_PATTERN) {
2992: MatZeroEntries(aij->B);
2993: }
2994: }
2995: Bdisassembled = PETSC_FALSE;
2996: if (!aij->B) {
2997: MatCreate(PETSC_COMM_SELF,&aij->B);
2998: PetscLogObjectParent((PetscObject)C,(PetscObject)aij->B);
2999: MatSetSizes(aij->B,C->rmap->n,C->cmap->N,C->rmap->n,C->cmap->N);
3000: MatSetBlockSizesFromMats(aij->B,B,B);
3001: MatSetType(aij->B,MATSEQAIJ);
3002: Bdisassembled = PETSC_TRUE;
3003: }
3004: if (B) {
3005: Baij = (Mat_SeqAIJ*)B->data;
3006: if (pattern == DIFFERENT_NONZERO_PATTERN) {
3007: PetscMalloc1(B->rmap->n,&nz);
3008: for (i=0; i<B->rmap->n; i++) {
3009: nz[i] = Baij->i[i+1] - Baij->i[i];
3010: }
3011: MatSeqAIJSetPreallocation(aij->B,0,nz);
3012: PetscFree(nz);
3013: }
3015: PetscLayoutGetRange(C->rmap,&rstart,&rend);
3016: shift = rend-rstart;
3017: count = 0;
3018: rowindices = NULL;
3019: colindices = NULL;
3020: if (rowemb) {
3021: ISGetIndices(rowemb,&rowindices);
3022: }
3023: if (ocolemb) {
3024: ISGetIndices(ocolemb,&colindices);
3025: }
3026: for (i=0; i<B->rmap->n; i++) {
3027: PetscInt row;
3028: row = i;
3029: if (rowindices) row = rowindices[i];
3030: for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
3031: col = Baij->j[count];
3032: if (colindices) col = colindices[col];
3033: if (Bdisassembled && col>=rstart) col += shift;
3034: v = Baij->a[count];
3035: MatSetValues(aij->B,1,&row,1,&col,&v,INSERT_VALUES);
3036: ++count;
3037: }
3038: }
3039: /* No assembly for aij->B is necessary. */
3040: /* FIXME: set aij->B's nonzerostate correctly. */
3041: } else {
3042: MatSetUp(aij->B);
3043: }
3044: C->preallocated = PETSC_TRUE;
3045: C->was_assembled = PETSC_FALSE;
3046: C->assembled = PETSC_FALSE;
3047: /*
3048: C will need to be assembled so that aij->B can be compressed into local form in MatSetUpMultiply_MPIAIJ().
3049: Furthermore, its nonzerostate will need to be based on that of aij->A's and aij->B's.
3050: */
3051: return(0);
3052: }
3054: /*
3055: B uses local indices with column indices ranging between 0 and N-n; they must be interpreted using garray.
3056: */
3057: PetscErrorCode MatGetSeqMats_MPIAIJ(Mat C,Mat *A,Mat *B)
3058: {
3059: Mat_MPIAIJ *aij = (Mat_MPIAIJ*)C->data;
3064: /* FIXME: make sure C is assembled */
3065: *A = aij->A;
3066: *B = aij->B;
3067: /* Note that we don't incref *A and *B, so be careful! */
3068: return(0);
3069: }
3071: /*
3072: Extract MPI submatrices encoded by pairs of IS that may live on subcomms of C.
3073: NOT SCALABLE due to the use of ISGetNonlocalIS() (see below).
3074: */
3075: PetscErrorCode MatCreateSubMatricesMPI_MPIXAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[],
3076: PetscErrorCode(*getsubmats_seq)(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat**),
3077: PetscErrorCode(*getlocalmats)(Mat,Mat*,Mat*),
3078: PetscErrorCode(*setseqmat)(Mat,IS,IS,MatStructure,Mat),
3079: PetscErrorCode(*setseqmats)(Mat,IS,IS,IS,MatStructure,Mat,Mat))
3080: {
3082: PetscMPIInt isize,flag;
3083: PetscInt i,ii,cismax,ispar;
3084: Mat *A,*B;
3085: IS *isrow_p,*iscol_p,*cisrow,*ciscol,*ciscol_p;
3088: if (!ismax) return(0);
3090: for (i = 0, cismax = 0; i < ismax; ++i) {
3091: PetscMPIInt isize;
3092: MPI_Comm_compare(((PetscObject)isrow[i])->comm,((PetscObject)iscol[i])->comm,&flag);
3093: if (flag != MPI_IDENT) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Row and column index sets must have the same communicator");
3094: MPI_Comm_size(((PetscObject)isrow[i])->comm, &isize);
3095: if (isize > 1) ++cismax;
3096: }
3098: /*
3099: If cismax is zero on all C's ranks, then and only then can we use purely sequential matrix extraction.
3100: ispar counts the number of parallel ISs across C's comm.
3101: */
3102: MPIU_Allreduce(&cismax,&ispar,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));
3103: if (!ispar) { /* Sequential ISs only across C's comm, so can call the sequential matrix extraction subroutine. */
3104: (*getsubmats_seq)(C,ismax,isrow,iscol,scall,submat);
3105: return(0);
3106: }
3108: /* if (ispar) */
3109: /*
3110: Construct the "complements" -- the off-processor indices -- of the iscol ISs for parallel ISs only.
3111: These are used to extract the off-diag portion of the resulting parallel matrix.
3112: The row IS for the off-diag portion is the same as for the diag portion,
3113: so we merely alias (without increfing) the row IS, while skipping those that are sequential.
3114: */
3115: PetscMalloc2(cismax,&cisrow,cismax,&ciscol);
3116: PetscMalloc1(cismax,&ciscol_p);
3117: for (i = 0, ii = 0; i < ismax; ++i) {
3118: MPI_Comm_size(((PetscObject)isrow[i])->comm,&isize);
3119: if (isize > 1) {
3120: /*
3121: TODO: This is the part that's ***NOT SCALABLE***.
3122: To fix this we need to extract just the indices of C's nonzero columns
3123: that lie on the intersection of isrow[i] and ciscol[ii] -- the nonlocal
3124: part of iscol[i] -- without actually computing ciscol[ii]. This also has
3125: to be done without serializing on the IS list, so, most likely, it is best
3126: done by rewriting MatCreateSubMatrices_MPIAIJ() directly.
3127: */
3128: ISGetNonlocalIS(iscol[i],&(ciscol[ii]));
3129: /* Now we have to
3130: (a) make sure ciscol[ii] is sorted, since, even if the off-proc indices
3131: were sorted on each rank, concatenated they might no longer be sorted;
3132: (b) Use ISSortPermutation() to construct ciscol_p, the mapping from the
3133: indices in the nondecreasing order to the original index positions.
3134: If ciscol[ii] is strictly increasing, the permutation IS is NULL.
3135: */
3136: ISSortPermutation(ciscol[ii],PETSC_FALSE,ciscol_p+ii);
3137: ISSort(ciscol[ii]);
3138: ++ii;
3139: }
3140: }
3141: PetscMalloc2(ismax,&isrow_p,ismax,&iscol_p);
3142: for (i = 0, ii = 0; i < ismax; ++i) {
3143: PetscInt j,issize;
3144: const PetscInt *indices;
3146: /*
3147: Permute the indices into a nondecreasing order. Reject row and col indices with duplicates.
3148: */
3149: ISSortPermutation(isrow[i],PETSC_FALSE,isrow_p+i);
3150: ISSort(isrow[i]);
3151: ISGetLocalSize(isrow[i],&issize);
3152: ISGetIndices(isrow[i],&indices);
3153: for (j = 1; j < issize; ++j) {
3154: if (indices[j] == indices[j-1]) {
3155: 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]);
3156: }
3157: }
3158: ISRestoreIndices(isrow[i],&indices);
3161: ISSortPermutation(iscol[i],PETSC_FALSE,iscol_p+i);
3162: ISSort(iscol[i]);
3163: ISGetLocalSize(iscol[i],&issize);
3164: ISGetIndices(iscol[i],&indices);
3165: for (j = 1; j < issize; ++j) {
3166: if (indices[j-1] == indices[j]) {
3167: 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]);
3168: }
3169: }
3170: ISRestoreIndices(iscol[i],&indices);
3171: MPI_Comm_size(((PetscObject)isrow[i])->comm,&isize);
3172: if (isize > 1) {
3173: cisrow[ii] = isrow[i];
3174: ++ii;
3175: }
3176: }
3177: /*
3178: Allocate the necessary arrays to hold the resulting parallel matrices as well as the intermediate
3179: array of sequential matrices underlying the resulting parallel matrices.
3180: Which arrays to allocate is based on the value of MatReuse scall and whether ISs are sorted and/or
3181: contain duplicates.
3183: There are as many diag matrices as there are original index sets. There are only as many parallel
3184: and off-diag matrices, as there are parallel (comm size > 1) index sets.
3186: ARRAYS that can hold Seq matrices get allocated in any event -- either here or by getsubmats_seq():
3187: - If the array of MPI matrices already exists and is being reused, we need to allocate the array
3188: and extract the underlying seq matrices into it to serve as placeholders, into which getsubmats_seq
3189: will deposite the extracted diag and off-diag parts. Thus, we allocate the A&B arrays and fill them
3190: with A[i] and B[ii] extracted from the corresponding MPI submat.
3191: - However, if the rows, A's column indices or B's column indices are not sorted, the extracted A[i] & B[ii]
3192: will have a different order from what getsubmats_seq expects. To handle this case -- indicated
3193: by a nonzero isrow_p[i], iscol_p[i], or ciscol_p[ii] -- we duplicate A[i] --> AA[i], B[ii] --> BB[ii]
3194: (retrieve composed AA[i] or BB[ii]) and reuse them here. AA[i] and BB[ii] are then used to permute its
3195: values into A[i] and B[ii] sitting inside the corresponding submat.
3196: - If no reuse is taking place then getsubmats_seq will allocate the A&B arrays and create the corresponding
3197: A[i], B[ii], AA[i] or BB[ii] matrices.
3198: */
3199: /* Parallel matrix array is allocated here only if no reuse is taking place. If reused, it is passed in by the caller. */
3200: if (scall == MAT_INITIAL_MATRIX) {
3201: PetscMalloc1(ismax,submat);
3202: }
3204: /* Now obtain the sequential A and B submatrices separately. */
3205: /* scall=MAT_REUSE_MATRIX is not handled yet, because getsubmats_seq() requires reuse of A and B */
3206: (*getsubmats_seq)(C,ismax,isrow,iscol,MAT_INITIAL_MATRIX,&A);
3207: (*getsubmats_seq)(C,cismax,cisrow,ciscol,MAT_INITIAL_MATRIX,&B);
3209: /*
3210: If scall == MAT_REUSE_MATRIX AND the permutations are NULL, we are done, since the sequential
3211: matrices A & B have been extracted directly into the parallel matrices containing them, or
3212: simply into the sequential matrix identical with the corresponding A (if isize == 1).
3213: Note that in that case colmap doesn't need to be rebuilt, since the matrices are expected
3214: to have the same sparsity pattern.
3215: Otherwise, A and/or B have to be properly embedded into C's index spaces and the correct colmap
3216: must be constructed for C. This is done by setseqmat(s).
3217: */
3218: for (i = 0, ii = 0; i < ismax; ++i) {
3219: /*
3220: TODO: cache ciscol, permutation ISs and maybe cisrow? What about isrow & iscol?
3221: That way we can avoid sorting and computing permutations when reusing.
3222: To this end:
3223: - remove the old cache, if it exists, when extracting submatrices with MAT_INITIAL_MATRIX
3224: - if caching arrays to hold the ISs, make and compose a container for them so that it can
3225: be destroyed upon destruction of C (use PetscContainerUserDestroy() to clear out the contents).
3226: */
3227: MatStructure pattern;
3228: pattern = DIFFERENT_NONZERO_PATTERN;
3230: MPI_Comm_size(((PetscObject)isrow[i])->comm,&isize);
3231: /* Construct submat[i] from the Seq pieces A (and B, if necessary). */
3232: if (isize > 1) {
3233: if (scall == MAT_INITIAL_MATRIX) {
3234: MatCreate(((PetscObject)isrow[i])->comm,(*submat)+i);
3235: MatSetSizes((*submat)[i],A[i]->rmap->n,A[i]->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
3236: MatSetType((*submat)[i],MATMPIAIJ);
3237: PetscLayoutSetUp((*submat)[i]->rmap);
3238: PetscLayoutSetUp((*submat)[i]->cmap);
3239: }
3240: /*
3241: For each parallel isrow[i], insert the extracted sequential matrices into the parallel matrix.
3242: */
3243: {
3244: Mat AA,BB;
3245: AA = A[i];
3246: BB = B[ii];
3247: if (AA || BB) {
3248: setseqmats((*submat)[i],isrow_p[i],iscol_p[i],ciscol_p[ii],pattern,AA,BB);
3249: MatAssemblyBegin((*submat)[i],MAT_FINAL_ASSEMBLY);
3250: MatAssemblyEnd((*submat)[i],MAT_FINAL_ASSEMBLY);
3251: }
3253: MatDestroy(&AA);
3254: }
3255: ISDestroy(ciscol+ii);
3256: ISDestroy(ciscol_p+ii);
3257: ++ii;
3258: } else { /* if (isize == 1) */
3259: if (scall == MAT_REUSE_MATRIX) {
3260: MatDestroy(&(*submat)[i]);
3261: }
3262: if (isrow_p[i] || iscol_p[i]) {
3263: MatDuplicate(A[i],MAT_DO_NOT_COPY_VALUES,(*submat)+i);
3264: setseqmat((*submat)[i],isrow_p[i],iscol_p[i],pattern,A[i]);
3265: /* Otherwise A is extracted straight into (*submats)[i]. */
3266: /* TODO: Compose A[i] on (*submat([i] for future use, if ((isrow_p[i] || iscol_p[i]) && MAT_INITIAL_MATRIX). */
3267: MatDestroy(A+i);
3268: } else (*submat)[i] = A[i];
3269: }
3270: ISDestroy(&isrow_p[i]);
3271: ISDestroy(&iscol_p[i]);
3272: }
3273: PetscFree2(cisrow,ciscol);
3274: PetscFree2(isrow_p,iscol_p);
3275: PetscFree(ciscol_p);
3276: PetscFree(A);
3277: MatDestroySubMatrices(cismax,&B);
3278: return(0);
3279: }
3281: PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
3282: {
3286: MatCreateSubMatricesMPI_MPIXAIJ(C,ismax,isrow,iscol,scall,submat,MatCreateSubMatrices_MPIAIJ,MatGetSeqMats_MPIAIJ,MatSetSeqMat_SeqAIJ,MatSetSeqMats_MPIAIJ);
3287: return(0);
3288: }