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