xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 0cd88d33dca7e1f18a10cbb6fcb08f83d068c5f4)
1 
2 /*
3   Defines projective product routines where A is a MPIAIJ matrix
4           C = P^T * A * P
5 */
6 
7 #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8 #include <../src/mat/utils/freespace.h>
9 #include <../src/mat/impls/aij/mpi/mpiaij.h>
10 #include <petscbt.h>
11 #include <petsctime.h>
12 
13 /* #define PTAP_PROFILE */
14 
15 PetscErrorCode MatView_MPIAIJ_PtAP(Mat A,PetscViewer viewer)
16 {
17   PetscErrorCode    ierr;
18   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data;
19   Mat_PtAPMPI       *ptap=a->ptap;
20   PetscBool         iascii;
21   PetscViewerFormat format;
22 
23   PetscFunctionBegin;
24   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
25   if (iascii) {
26     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
27     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
28       if (ptap->algType == 0) {
29         ierr = PetscViewerASCIIPrintf(viewer,"using scalable MatPtAP() implementation\n");CHKERRQ(ierr);
30       } else if (ptap->algType == 1) {
31         ierr = PetscViewerASCIIPrintf(viewer,"using nonscalable MatPtAP() implementation\n");CHKERRQ(ierr);
32       }
33     }
34   }
35   ierr = (ptap->view)(A,viewer);CHKERRQ(ierr);
36   PetscFunctionReturn(0);
37 }
38 
39 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
40 {
41   PetscErrorCode ierr;
42   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
43   Mat_PtAPMPI    *ptap=a->ptap;
44 
45   PetscFunctionBegin;
46   if (ptap) {
47     Mat_Merge_SeqsToMPI *merge=ptap->merge;
48     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
49     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
50     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
51     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
52     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
53     ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
54     ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
55     if (ptap->AP_loc) { /* used by alg_rap */
56       Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
57       ierr = PetscFree(ap->i);CHKERRQ(ierr);
58       ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
59       ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
60     } else { /* used by alg_ptap */
61       ierr = PetscFree(ptap->api);CHKERRQ(ierr);
62       ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
63     }
64     ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
65     ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
66     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
67 
68     if (merge) { /* used by alg_ptap */
69       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
70       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
71       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
72       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
73       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
74       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
75       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
76       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
77       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
78       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
79       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
80       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
81       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
82       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
83     }
84 
85     ierr = ptap->destroy(A);CHKERRQ(ierr);
86     ierr = PetscFree(ptap);CHKERRQ(ierr);
87   }
88   PetscFunctionReturn(0);
89 }
90 
91 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
92 {
93   PetscErrorCode ierr;
94   PetscBool      flg;
95   MPI_Comm       comm;
96 #if !defined(PETSC_HAVE_HYPRE)
97   const char          *algTypes[2] = {"scalable","nonscalable"};
98   PetscInt            nalg=2;
99 #else
100   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
101   PetscInt            nalg=3;
102 #endif
103   PetscInt            pN=P->cmap->N,alg=1; /* set default algorithm */
104 
105   PetscFunctionBegin;
106   /* check if matrix local sizes are compatible */
107   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
108   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
109   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
110 
111   if (scall == MAT_INITIAL_MATRIX) {
112     /* pick an algorithm */
113     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
114     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
115     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
116     ierr = PetscOptionsEnd();CHKERRQ(ierr);
117 
118     if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */
119       MatInfo     Ainfo,Pinfo;
120       PetscInt    nz_local;
121       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
122 
123       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
124       ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
125       nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
126 
127       if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
128       ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
129 
130       if (alg_scalable) {
131         alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
132       }
133     }
134 
135     switch (alg) {
136     case 1:
137       /* do R=P^T locally, then C=R*A*P */
138       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
139       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
140       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
141       break;
142 #if defined(PETSC_HAVE_HYPRE)
143     case 2:
144       /* Use boomerAMGBuildCoarseOperator */
145       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
146       PetscFunctionReturn(0);
147       break;
148 #endif
149     default:
150       /* do R=P^T locally, then C=R*A*P */
151       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
152       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr);
153       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
154       break;
155     }
156   }
157   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
158   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
159   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
160   PetscFunctionReturn(0);
161 }
162 
163 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C)
164 {
165   PetscErrorCode    ierr;
166   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
167   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
168   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
169   Mat_PtAPMPI       *ptap = c->ptap;
170   Mat               AP_loc,C_loc,C_oth;
171   PetscInt          i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz;
172   PetscScalar       *apa;
173   const PetscInt    *cols;
174   const PetscScalar *vals;
175 
176   PetscFunctionBegin;
177   ierr = MatZeroEntries(C);CHKERRQ(ierr);
178 
179   /* 1) get R = Pd^T,Ro = Po^T */
180   if (ptap->reuse == MAT_REUSE_MATRIX) {
181     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
182     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
183   }
184 
185   /* 2) get AP_loc */
186   AP_loc = ptap->AP_loc;
187   ap = (Mat_SeqAIJ*)AP_loc->data;
188 
189   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
190   /*-----------------------------------------------------*/
191   if (ptap->reuse == MAT_REUSE_MATRIX) {
192     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
193     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
194     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
195   }
196 
197   /* 2-2) compute numeric A_loc*P - dominating part */
198   /* ---------------------------------------------- */
199   /* get data from symbolic products */
200   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
201   if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
202   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ptap->P_oth is NULL. Cannot proceed!");
203 
204   api   = ap->i;
205   apj   = ap->j;
206   for (i=0; i<am; i++) {
207     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
208     apnz = api[i+1] - api[i];
209     apa = ap->a + api[i];
210     ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr);
211     AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa);
212     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
213   }
214 
215   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
216   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
217   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
218 
219   C_loc = ptap->C_loc;
220   C_oth = ptap->C_oth;
221 
222   /* add C_loc and Co to to C */
223   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
224 
225   /* C_loc -> C */
226   cm    = C_loc->rmap->N;
227   c_seq = (Mat_SeqAIJ*)C_loc->data;
228   cols = c_seq->j;
229   vals = c_seq->a;
230 
231   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
232   /* when there are no off-processor parts.  */
233   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
234   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
235   /* a table, and other, more complex stuff has to be done. */
236   if (C->assembled) {
237     C->was_assembled = PETSC_TRUE;
238     C->assembled     = PETSC_FALSE;
239   }
240   if (C->was_assembled) {
241     for (i=0; i<cm; i++) {
242       ncols = c_seq->i[i+1] - c_seq->i[i];
243       row = rstart + i;
244       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
245       cols += ncols; vals += ncols;
246     }
247   } else {
248     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
249   }
250 
251   /* Co -> C, off-processor part */
252   cm = C_oth->rmap->N;
253   c_seq = (Mat_SeqAIJ*)C_oth->data;
254   cols = c_seq->j;
255   vals = c_seq->a;
256   for (i=0; i<cm; i++) {
257     ncols = c_seq->i[i+1] - c_seq->i[i];
258     row = p->garray[i];
259     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
260     cols += ncols; vals += ncols;
261   }
262   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
263   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
264 
265   ptap->reuse = MAT_REUSE_MATRIX;
266   PetscFunctionReturn(0);
267 }
268 
269 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
270 {
271   PetscErrorCode      ierr;
272   Mat_PtAPMPI         *ptap;
273   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
274   MPI_Comm            comm;
275   PetscMPIInt         size,rank;
276   Mat                 Cmpi,P_loc,P_oth;
277   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
278   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
279   PetscInt            *lnk,i,k,pnz,row,nsend;
280   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
281   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
282   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
283   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
284   MPI_Request         *swaits,*rwaits;
285   MPI_Status          *sstatus,rstatus;
286   PetscLayout         rowmap;
287   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
288   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
289   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi;
290   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
291   PetscScalar         *apv;
292   PetscTable          ta;
293   MatType             mtype;
294 #if defined(PETSC_USE_INFO)
295   PetscReal           apfill;
296 #endif
297 
298   PetscFunctionBegin;
299   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
300   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
301   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
302 
303   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
304 
305   /* create symbolic parallel matrix Cmpi */
306   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
307   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
308   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
309 
310   /* create struct Mat_PtAPMPI and attached it to C later */
311   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
312   ptap->reuse = MAT_INITIAL_MATRIX;
313   ptap->algType = 0;
314 
315   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
316   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
317   /* get P_loc by taking all local rows of P */
318   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
319 
320   ptap->P_loc = P_loc;
321   ptap->P_oth = P_oth;
322 
323   /* (0) compute Rd = Pd^T, Ro = Po^T  */
324   /* --------------------------------- */
325   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
326   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
327 
328   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
329   /* ----------------------------------------------------------------- */
330   p_loc  = (Mat_SeqAIJ*)P_loc->data;
331   if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data;
332 
333   /* create and initialize a linked list */
334   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
335   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
336   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
337   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
338 
339   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
340 
341   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
342   if (ao) {
343     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
344   } else {
345     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
346   }
347   current_space = free_space;
348   nspacedouble  = 0;
349 
350   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
351   api[0] = 0;
352   for (i=0; i<am; i++) {
353     /* diagonal portion: Ad[i,:]*P */
354     ai = ad->i; pi = p_loc->i;
355     nzi = ai[i+1] - ai[i];
356     aj  = ad->j + ai[i];
357     for (j=0; j<nzi; j++) {
358       row  = aj[j];
359       pnz  = pi[row+1] - pi[row];
360       Jptr = p_loc->j + pi[row];
361       /* add non-zero cols of P into the sorted linked list lnk */
362       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
363     }
364     /* off-diagonal portion: Ao[i,:]*P */
365     if (ao) {
366       ai = ao->i; pi = p_oth->i;
367       nzi = ai[i+1] - ai[i];
368       aj  = ao->j + ai[i];
369       for (j=0; j<nzi; j++) {
370         row  = aj[j];
371         pnz  = pi[row+1] - pi[row];
372         Jptr = p_oth->j + pi[row];
373         ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
374       }
375     }
376     apnz     = lnk[0];
377     api[i+1] = api[i] + apnz;
378 
379     /* if free space is not available, double the total space in the list */
380     if (current_space->local_remaining<apnz) {
381       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
382       nspacedouble++;
383     }
384 
385     /* Copy data into free space, then initialize lnk */
386     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
387 
388     current_space->array           += apnz;
389     current_space->local_used      += apnz;
390     current_space->local_remaining -= apnz;
391   }
392   /* Allocate space for apj and apv, initialize apj, and */
393   /* destroy list of free space and other temporary array(s) */
394   ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
395   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
396   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
397 
398   /* Create AP_loc for reuse */
399   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
400 
401 #if defined(PETSC_USE_INFO)
402   if (ao) {
403     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
404   } else {
405     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
406   }
407   ptap->AP_loc->info.mallocs           = nspacedouble;
408   ptap->AP_loc->info.fill_ratio_given  = fill;
409   ptap->AP_loc->info.fill_ratio_needed = apfill;
410 
411   if (api[am]) {
412     ierr = PetscInfo3(ptap->AP_loc,"Scalable algorithm, AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
413     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
414   } else {
415     ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
416   }
417 #endif
418 
419   /* (2-1) compute symbolic Co = Ro*AP_loc  */
420   /* ------------------------------------ */
421   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
422 
423   /* (3) send coj of C_oth to other processors  */
424   /* ------------------------------------------ */
425   /* determine row ownership */
426   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
427   rowmap->n  = pn;
428   rowmap->bs = 1;
429   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
430   owners = rowmap->range;
431 
432   /* determine the number of messages to send, their lengths */
433   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
434   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
435   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
436 
437   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
438   coi   = c_oth->i; coj = c_oth->j;
439   con   = ptap->C_oth->rmap->n;
440   proc  = 0;
441   for (i=0; i<con; i++) {
442     while (prmap[i] >= owners[proc+1]) proc++;
443     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
444     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
445   }
446 
447   len          = 0; /* max length of buf_si[], see (4) */
448   owners_co[0] = 0;
449   nsend        = 0;
450   for (proc=0; proc<size; proc++) {
451     owners_co[proc+1] = owners_co[proc] + len_si[proc];
452     if (len_s[proc]) {
453       nsend++;
454       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
455       len         += len_si[proc];
456     }
457   }
458 
459   /* determine the number and length of messages to receive for coi and coj  */
460   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
461   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
462 
463   /* post the Irecv and Isend of coj */
464   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
465   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
466   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
467   for (proc=0, k=0; proc<size; proc++) {
468     if (!len_s[proc]) continue;
469     i    = owners_co[proc];
470     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
471     k++;
472   }
473 
474   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
475   /* ---------------------------------------- */
476   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
477   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
478 
479   /* receives coj are complete */
480   for (i=0; i<nrecv; i++) {
481     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
482   }
483   ierr = PetscFree(rwaits);CHKERRQ(ierr);
484   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
485 
486   /* add received column indices into ta to update Crmax */
487   for (k=0; k<nrecv; k++) {/* k-th received message */
488     Jptr = buf_rj[k];
489     for (j=0; j<len_r[k]; j++) {
490       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
491     }
492   }
493   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
494   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
495 
496   /* (4) send and recv coi */
497   /*-----------------------*/
498   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
499   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
500   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
501   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
502   for (proc=0,k=0; proc<size; proc++) {
503     if (!len_s[proc]) continue;
504     /* form outgoing message for i-structure:
505          buf_si[0]:                 nrows to be sent
506                [1:nrows]:           row index (global)
507                [nrows+1:2*nrows+1]: i-structure index
508     */
509     /*-------------------------------------------*/
510     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
511     buf_si_i    = buf_si + nrows+1;
512     buf_si[0]   = nrows;
513     buf_si_i[0] = 0;
514     nrows       = 0;
515     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
516       nzi = coi[i+1] - coi[i];
517       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
518       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
519       nrows++;
520     }
521     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
522     k++;
523     buf_si += len_si[proc];
524   }
525   for (i=0; i<nrecv; i++) {
526     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
527   }
528   ierr = PetscFree(rwaits);CHKERRQ(ierr);
529   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
530 
531   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
532   ierr = PetscFree(len_ri);CHKERRQ(ierr);
533   ierr = PetscFree(swaits);CHKERRQ(ierr);
534   ierr = PetscFree(buf_s);CHKERRQ(ierr);
535 
536   /* (5) compute the local portion of Cmpi      */
537   /* ------------------------------------------ */
538   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
539   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
540   current_space = free_space;
541 
542   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
543   for (k=0; k<nrecv; k++) {
544     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
545     nrows       = *buf_ri_k[k];
546     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
547     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
548   }
549 
550   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
551   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
552   for (i=0; i<pn; i++) {
553     /* add C_loc into Cmpi */
554     nzi  = c_loc->i[i+1] - c_loc->i[i];
555     Jptr = c_loc->j + c_loc->i[i];
556     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
557 
558     /* add received col data into lnk */
559     for (k=0; k<nrecv; k++) { /* k-th received message */
560       if (i == *nextrow[k]) { /* i-th row */
561         nzi  = *(nextci[k]+1) - *nextci[k];
562         Jptr = buf_rj[k] + *nextci[k];
563         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
564         nextrow[k]++; nextci[k]++;
565       }
566     }
567     nzi = lnk[0];
568 
569     /* copy data into free space, then initialize lnk */
570     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
571     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
572   }
573   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
574   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
575   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
576 
577   /* local sizes and preallocation */
578   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
579   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
580   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
581   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
582 
583   /* members in merge */
584   ierr = PetscFree(id_r);CHKERRQ(ierr);
585   ierr = PetscFree(len_r);CHKERRQ(ierr);
586   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
587   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
588   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
589   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
590   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
591 
592   /* attach the supporting struct to Cmpi for reuse */
593   c = (Mat_MPIAIJ*)Cmpi->data;
594   c->ptap         = ptap;
595   ptap->duplicate = Cmpi->ops->duplicate;
596   ptap->destroy   = Cmpi->ops->destroy;
597   ptap->view      = Cmpi->ops->view;
598 
599   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
600   Cmpi->assembled        = PETSC_FALSE;
601   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
602   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
603   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
604   *C                     = Cmpi;
605   PetscFunctionReturn(0);
606 }
607 
608 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
609 {
610   PetscErrorCode      ierr;
611   Mat_PtAPMPI         *ptap;
612   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
613   MPI_Comm            comm;
614   PetscMPIInt         size,rank;
615   Mat                 Cmpi;
616   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
617   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
618   PetscInt            *lnk,i,k,pnz,row,nsend;
619   PetscBT             lnkbt;
620   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
621   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
622   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
623   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
624   MPI_Request         *swaits,*rwaits;
625   MPI_Status          *sstatus,rstatus;
626   PetscLayout         rowmap;
627   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
628   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
629   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
630   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
631   PetscScalar         *apv;
632   PetscTable          ta;
633   MatType             mtype;
634 #if defined(PETSC_USE_INFO)
635   PetscReal           apfill;
636 #endif
637 
638   PetscFunctionBegin;
639   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
640   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
641   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
642 
643   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
644 
645   /* create symbolic parallel matrix Cmpi */
646   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
647   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
648   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
649 
650   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
651   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
652 
653   /* create struct Mat_PtAPMPI and attached it to C later */
654   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
655   ptap->reuse = MAT_INITIAL_MATRIX;
656   ptap->algType = 1;
657 
658   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
659   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
660   /* get P_loc by taking all local rows of P */
661   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
662 
663   /* (0) compute Rd = Pd^T, Ro = Po^T  */
664   /* --------------------------------- */
665   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
666   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
667 
668   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
669   /* ----------------------------------------------------------------- */
670   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
671   if (ptap->P_oth) p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
672 
673   /* create and initialize a linked list */
674   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
675   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
676   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
677   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
678   /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */
679 
680   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
681 
682   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
683   if (ao) {
684     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
685   } else {
686     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
687   }
688   current_space = free_space;
689   nspacedouble  = 0;
690 
691   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
692   api[0] = 0;
693   for (i=0; i<am; i++) {
694     /* diagonal portion: Ad[i,:]*P */
695     ai = ad->i; pi = p_loc->i;
696     nzi = ai[i+1] - ai[i];
697     aj  = ad->j + ai[i];
698     for (j=0; j<nzi; j++) {
699       row  = aj[j];
700       pnz  = pi[row+1] - pi[row];
701       Jptr = p_loc->j + pi[row];
702       /* add non-zero cols of P into the sorted linked list lnk */
703       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
704     }
705     /* off-diagonal portion: Ao[i,:]*P */
706     if (ao) {
707       ai = ao->i; pi = p_oth->i;
708       nzi = ai[i+1] - ai[i];
709       aj  = ao->j + ai[i];
710       for (j=0; j<nzi; j++) {
711         row  = aj[j];
712         pnz  = pi[row+1] - pi[row];
713         Jptr = p_oth->j + pi[row];
714         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
715       }
716     }
717     apnz     = lnk[0];
718     api[i+1] = api[i] + apnz;
719     if (ap_rmax < apnz) ap_rmax = apnz;
720 
721     /* if free space is not available, double the total space in the list */
722     if (current_space->local_remaining<apnz) {
723       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
724       nspacedouble++;
725     }
726 
727     /* Copy data into free space, then initialize lnk */
728     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
729 
730     current_space->array           += apnz;
731     current_space->local_used      += apnz;
732     current_space->local_remaining -= apnz;
733   }
734   /* Allocate space for apj and apv, initialize apj, and */
735   /* destroy list of free space and other temporary array(s) */
736   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
737   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
738   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
739 
740   /* Create AP_loc for reuse */
741   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
742 
743 #if defined(PETSC_USE_INFO)
744   if (ao) {
745     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
746   } else {
747     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
748   }
749   ptap->AP_loc->info.mallocs           = nspacedouble;
750   ptap->AP_loc->info.fill_ratio_given  = fill;
751   ptap->AP_loc->info.fill_ratio_needed = apfill;
752 
753   if (api[am]) {
754     ierr = PetscInfo3(ptap->AP_loc,"Nonscalable algorithm, AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
755     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
756   } else {
757     ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
758   }
759 #endif
760 
761   /* (2-1) compute symbolic Co = Ro*AP_loc  */
762   /* ------------------------------------ */
763   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
764 
765   /* (3) send coj of C_oth to other processors  */
766   /* ------------------------------------------ */
767   /* determine row ownership */
768   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
769   rowmap->n  = pn;
770   rowmap->bs = 1;
771   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
772   owners = rowmap->range;
773 
774   /* determine the number of messages to send, their lengths */
775   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
776   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
777   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
778 
779   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
780   coi   = c_oth->i; coj = c_oth->j;
781   con   = ptap->C_oth->rmap->n;
782   proc  = 0;
783   for (i=0; i<con; i++) {
784     while (prmap[i] >= owners[proc+1]) proc++;
785     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
786     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
787   }
788 
789   len          = 0; /* max length of buf_si[], see (4) */
790   owners_co[0] = 0;
791   nsend        = 0;
792   for (proc=0; proc<size; proc++) {
793     owners_co[proc+1] = owners_co[proc] + len_si[proc];
794     if (len_s[proc]) {
795       nsend++;
796       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
797       len         += len_si[proc];
798     }
799   }
800 
801   /* determine the number and length of messages to receive for coi and coj  */
802   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
803   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
804 
805   /* post the Irecv and Isend of coj */
806   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
807   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
808   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
809   for (proc=0, k=0; proc<size; proc++) {
810     if (!len_s[proc]) continue;
811     i    = owners_co[proc];
812     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
813     k++;
814   }
815 
816   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
817   /* ---------------------------------------- */
818   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
819   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
820 
821   /* receives coj are complete */
822   for (i=0; i<nrecv; i++) {
823     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
824   }
825   ierr = PetscFree(rwaits);CHKERRQ(ierr);
826   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
827 
828   /* add received column indices into ta to update Crmax */
829   for (k=0; k<nrecv; k++) {/* k-th received message */
830     Jptr = buf_rj[k];
831     for (j=0; j<len_r[k]; j++) {
832       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
833     }
834   }
835   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
836   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
837 
838   /* (4) send and recv coi */
839   /*-----------------------*/
840   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
841   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
842   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
843   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
844   for (proc=0,k=0; proc<size; proc++) {
845     if (!len_s[proc]) continue;
846     /* form outgoing message for i-structure:
847          buf_si[0]:                 nrows to be sent
848                [1:nrows]:           row index (global)
849                [nrows+1:2*nrows+1]: i-structure index
850     */
851     /*-------------------------------------------*/
852     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
853     buf_si_i    = buf_si + nrows+1;
854     buf_si[0]   = nrows;
855     buf_si_i[0] = 0;
856     nrows       = 0;
857     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
858       nzi = coi[i+1] - coi[i];
859       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
860       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
861       nrows++;
862     }
863     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
864     k++;
865     buf_si += len_si[proc];
866   }
867   for (i=0; i<nrecv; i++) {
868     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
869   }
870   ierr = PetscFree(rwaits);CHKERRQ(ierr);
871   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
872 
873   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
874   ierr = PetscFree(len_ri);CHKERRQ(ierr);
875   ierr = PetscFree(swaits);CHKERRQ(ierr);
876   ierr = PetscFree(buf_s);CHKERRQ(ierr);
877 
878   /* (5) compute the local portion of Cmpi      */
879   /* ------------------------------------------ */
880   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
881   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
882   current_space = free_space;
883 
884   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
885   for (k=0; k<nrecv; k++) {
886     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
887     nrows       = *buf_ri_k[k];
888     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
889     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
890   }
891 
892   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
893   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
894   for (i=0; i<pn; i++) {
895     /* add C_loc into Cmpi */
896     nzi  = c_loc->i[i+1] - c_loc->i[i];
897     Jptr = c_loc->j + c_loc->i[i];
898     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
899 
900     /* add received col data into lnk */
901     for (k=0; k<nrecv; k++) { /* k-th received message */
902       if (i == *nextrow[k]) { /* i-th row */
903         nzi  = *(nextci[k]+1) - *nextci[k];
904         Jptr = buf_rj[k] + *nextci[k];
905         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
906         nextrow[k]++; nextci[k]++;
907       }
908     }
909     nzi = lnk[0];
910 
911     /* copy data into free space, then initialize lnk */
912     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
913     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
914   }
915   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
916   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
917   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
918 
919   /* local sizes and preallocation */
920   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
921   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
922   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
923   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
924 
925   /* members in merge */
926   ierr = PetscFree(id_r);CHKERRQ(ierr);
927   ierr = PetscFree(len_r);CHKERRQ(ierr);
928   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
929   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
930   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
931   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
932   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
933 
934   /* attach the supporting struct to Cmpi for reuse */
935   c = (Mat_MPIAIJ*)Cmpi->data;
936   c->ptap         = ptap;
937   ptap->duplicate = Cmpi->ops->duplicate;
938   ptap->destroy   = Cmpi->ops->destroy;
939   ptap->view      = Cmpi->ops->view;
940   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
941 
942   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
943   Cmpi->assembled        = PETSC_FALSE;
944   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
945   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
946   *C                     = Cmpi;
947   PetscFunctionReturn(0);
948 }
949 
950 
951 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
952 {
953   PetscErrorCode    ierr;
954   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
955   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
956   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
957   Mat_PtAPMPI       *ptap = c->ptap;
958   Mat               AP_loc,C_loc,C_oth;
959   PetscInt          i,rstart,rend,cm,ncols,row;
960   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
961   PetscScalar       *apa;
962   const PetscInt    *cols;
963   const PetscScalar *vals;
964 
965   PetscFunctionBegin;
966   ierr = MatZeroEntries(C);CHKERRQ(ierr);
967   /* 1) get R = Pd^T,Ro = Po^T */
968   if (ptap->reuse == MAT_REUSE_MATRIX) {
969     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
970     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
971   }
972 
973   /* 2) get AP_loc */
974   AP_loc = ptap->AP_loc;
975   ap = (Mat_SeqAIJ*)AP_loc->data;
976 
977   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
978   /*-----------------------------------------------------*/
979   if (ptap->reuse == MAT_REUSE_MATRIX) {
980     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
981     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
982     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
983   }
984 
985   /* 2-2) compute numeric A_loc*P - dominating part */
986   /* ---------------------------------------------- */
987   /* get data from symbolic products */
988   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
989   if (ptap->P_oth) {
990     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
991   }
992   apa   = ptap->apa;
993   api   = ap->i;
994   apj   = ap->j;
995   for (i=0; i<am; i++) {
996     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
997     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
998     apnz = api[i+1] - api[i];
999     for (j=0; j<apnz; j++) {
1000       col = apj[j+api[i]];
1001       ap->a[j+ap->i[i]] = apa[col];
1002       apa[col] = 0.0;
1003     }
1004     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1005   }
1006 
1007   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
1008   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
1009   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
1010   C_loc = ptap->C_loc;
1011   C_oth = ptap->C_oth;
1012 
1013   /* add C_loc and Co to to C */
1014   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
1015 
1016   /* C_loc -> C */
1017   cm    = C_loc->rmap->N;
1018   c_seq = (Mat_SeqAIJ*)C_loc->data;
1019   cols = c_seq->j;
1020   vals = c_seq->a;
1021 
1022 
1023   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
1024   /* when there are no off-processor parts.  */
1025   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
1026   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
1027   /* a table, and other, more complex stuff has to be done. */
1028   if (C->assembled) {
1029     C->was_assembled = PETSC_TRUE;
1030     C->assembled     = PETSC_FALSE;
1031   }
1032   if (C->was_assembled) {
1033     for (i=0; i<cm; i++) {
1034       ncols = c_seq->i[i+1] - c_seq->i[i];
1035       row = rstart + i;
1036       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1037       cols += ncols; vals += ncols;
1038     }
1039   } else {
1040     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
1041   }
1042 
1043   /* Co -> C, off-processor part */
1044   cm = C_oth->rmap->N;
1045   c_seq = (Mat_SeqAIJ*)C_oth->data;
1046   cols = c_seq->j;
1047   vals = c_seq->a;
1048   for (i=0; i<cm; i++) {
1049     ncols = c_seq->i[i+1] - c_seq->i[i];
1050     row = p->garray[i];
1051     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1052     cols += ncols; vals += ncols;
1053   }
1054 
1055   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1056   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1057 
1058   ptap->reuse = MAT_REUSE_MATRIX;
1059   PetscFunctionReturn(0);
1060 }
1061