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