xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision df97dc6db963c80699cff554c068a8f1b2a115cb)
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   if (ptap->C_loc->ops->matmultnumeric) {
257     ierr = (*ptap->C_loc->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
258   } else {
259     ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
260   }
261   if (ptap->C_oth->ops->matmultnumeric) {
262     ierr = (*ptap->C_oth->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
263   } else {
264     ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
265   }
266 
267   C_loc = ptap->C_loc;
268   C_oth = ptap->C_oth;
269 
270   /* add C_loc and Co to to C */
271   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
272 
273   /* C_loc -> C */
274   cm    = C_loc->rmap->N;
275   c_seq = (Mat_SeqAIJ*)C_loc->data;
276   cols = c_seq->j;
277   vals = c_seq->a;
278   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_loc->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr);
279 
280   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
281   /* when there are no off-processor parts.  */
282   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
283   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
284   /* a table, and other, more complex stuff has to be done. */
285   if (C->assembled) {
286     C->was_assembled = PETSC_TRUE;
287     C->assembled     = PETSC_FALSE;
288   }
289   if (C->was_assembled) {
290     for (i=0; i<cm; i++) {
291       ncols = c_seq->i[i+1] - c_seq->i[i];
292       row = rstart + i;
293       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
294       cols += ncols; vals += ncols;
295     }
296   } else {
297     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
298   }
299   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_loc->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr);
300   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);
301 
302   /* Co -> C, off-processor part */
303   cm = C_oth->rmap->N;
304   c_seq = (Mat_SeqAIJ*)C_oth->data;
305   cols = c_seq->j;
306   vals = c_seq->a;
307   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_oth->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr);
308   for (i=0; i<cm; i++) {
309     ncols = c_seq->i[i+1] - c_seq->i[i];
310     row = p->garray[i];
311     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
312     cols += ncols; vals += ncols;
313   }
314   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
315   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
316 
317   ptap->reuse = MAT_REUSE_MATRIX;
318 
319   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_oth->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr);
320   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);
321 
322   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
323   if (ptap->freestruct) {
324     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
325   }
326   PetscFunctionReturn(0);
327 }
328 
329 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
330 {
331   PetscErrorCode      ierr;
332   Mat_APMPI           *ptap;
333   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
334   MPI_Comm            comm;
335   PetscMPIInt         size,rank;
336   Mat                 Cmpi,P_loc,P_oth;
337   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
338   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
339   PetscInt            *lnk,i,k,pnz,row,nsend;
340   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
341   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
342   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
343   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
344   MPI_Request         *swaits,*rwaits;
345   MPI_Status          *sstatus,rstatus;
346   PetscLayout         rowmap;
347   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
348   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
349   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi,nout;
350   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
351   PetscScalar         *apv;
352   PetscTable          ta;
353   MatType             mtype;
354 #if defined(PETSC_USE_INFO)
355   PetscReal           apfill;
356 #endif
357 
358   PetscFunctionBegin;
359   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
360   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
361   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
362 
363   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
364 
365   /* create symbolic parallel matrix Cmpi */
366   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
367   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
368   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
369 
370   /* create struct Mat_APMPI and attached it to C later */
371   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
372   ptap->reuse = MAT_INITIAL_MATRIX;
373   ptap->algType = 0;
374 
375   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
376   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
377   /* get P_loc by taking all local rows of P */
378   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
379 
380   ptap->P_loc = P_loc;
381   ptap->P_oth = P_oth;
382 
383   /* (0) compute Rd = Pd^T, Ro = Po^T  */
384   /* --------------------------------- */
385   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
386   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
387 
388   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
389   /* ----------------------------------------------------------------- */
390   p_loc  = (Mat_SeqAIJ*)P_loc->data;
391   if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data;
392 
393   /* create and initialize a linked list */
394   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
395   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
396   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
397   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
398 
399   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
400 
401   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
402   if (ao) {
403     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
404   } else {
405     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
406   }
407   current_space = free_space;
408   nspacedouble  = 0;
409 
410   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
411   api[0] = 0;
412   for (i=0; i<am; i++) {
413     /* diagonal portion: Ad[i,:]*P */
414     ai = ad->i; pi = p_loc->i;
415     nzi = ai[i+1] - ai[i];
416     aj  = ad->j + ai[i];
417     for (j=0; j<nzi; j++) {
418       row  = aj[j];
419       pnz  = pi[row+1] - pi[row];
420       Jptr = p_loc->j + pi[row];
421       /* add non-zero cols of P into the sorted linked list lnk */
422       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
423     }
424     /* off-diagonal portion: Ao[i,:]*P */
425     if (ao) {
426       ai = ao->i; pi = p_oth->i;
427       nzi = ai[i+1] - ai[i];
428       aj  = ao->j + ai[i];
429       for (j=0; j<nzi; j++) {
430         row  = aj[j];
431         pnz  = pi[row+1] - pi[row];
432         Jptr = p_oth->j + pi[row];
433         ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
434       }
435     }
436     apnz     = lnk[0];
437     api[i+1] = api[i] + apnz;
438 
439     /* if free space is not available, double the total space in the list */
440     if (current_space->local_remaining<apnz) {
441       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
442       nspacedouble++;
443     }
444 
445     /* Copy data into free space, then initialize lnk */
446     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
447 
448     current_space->array           += apnz;
449     current_space->local_used      += apnz;
450     current_space->local_remaining -= apnz;
451   }
452   /* Allocate space for apj and apv, initialize apj, and */
453   /* destroy list of free space and other temporary array(s) */
454   ierr = PetscCalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
455   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
456   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
457 
458   /* Create AP_loc for reuse */
459   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
460   ierr = MatSeqAIJCompactOutExtraColumns_SeqAIJ(ptap->AP_loc, &ptap->ltog);CHKERRQ(ierr);
461 
462 #if defined(PETSC_USE_INFO)
463   if (ao) {
464     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
465   } else {
466     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
467   }
468   ptap->AP_loc->info.mallocs           = nspacedouble;
469   ptap->AP_loc->info.fill_ratio_given  = fill;
470   ptap->AP_loc->info.fill_ratio_needed = apfill;
471 
472   if (api[am]) {
473     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);
474     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
475   } else {
476     ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
477   }
478 #endif
479 
480   /* (2-1) compute symbolic Co = Ro*AP_loc  */
481   /* ------------------------------------ */
482   ierr = MatSetOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr);
483   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
484 
485   /* (3) send coj of C_oth to other processors  */
486   /* ------------------------------------------ */
487   /* determine row ownership */
488   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
489   rowmap->n  = pn;
490   rowmap->bs = 1;
491   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
492   owners = rowmap->range;
493 
494   /* determine the number of messages to send, their lengths */
495   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
496   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
497   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
498 
499   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
500   coi   = c_oth->i; coj = c_oth->j;
501   con   = ptap->C_oth->rmap->n;
502   proc  = 0;
503   ierr = ISLocalToGlobalMappingApply(ptap->ltog,coi[con],coj,coj);CHKERRQ(ierr);
504   for (i=0; i<con; i++) {
505     while (prmap[i] >= owners[proc+1]) proc++;
506     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
507     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
508   }
509 
510   len          = 0; /* max length of buf_si[], see (4) */
511   owners_co[0] = 0;
512   nsend        = 0;
513   for (proc=0; proc<size; proc++) {
514     owners_co[proc+1] = owners_co[proc] + len_si[proc];
515     if (len_s[proc]) {
516       nsend++;
517       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
518       len         += len_si[proc];
519     }
520   }
521 
522   /* determine the number and length of messages to receive for coi and coj  */
523   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
524   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
525 
526   /* post the Irecv and Isend of coj */
527   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
528   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
529   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
530   for (proc=0, k=0; proc<size; proc++) {
531     if (!len_s[proc]) continue;
532     i    = owners_co[proc];
533     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
534     k++;
535   }
536 
537   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
538   /* ---------------------------------------- */
539   ierr = MatSetOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr);
540   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
541   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
542   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,c_loc->j);CHKERRQ(ierr);
543 
544   /* receives coj are complete */
545   for (i=0; i<nrecv; i++) {
546     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
547   }
548   ierr = PetscFree(rwaits);CHKERRQ(ierr);
549   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
550 
551   /* add received column indices into ta to update Crmax */
552   for (k=0; k<nrecv; k++) {/* k-th received message */
553     Jptr = buf_rj[k];
554     for (j=0; j<len_r[k]; j++) {
555       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
556     }
557   }
558   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
559   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
560 
561   /* (4) send and recv coi */
562   /*-----------------------*/
563   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
564   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
565   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
566   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
567   for (proc=0,k=0; proc<size; proc++) {
568     if (!len_s[proc]) continue;
569     /* form outgoing message for i-structure:
570          buf_si[0]:                 nrows to be sent
571                [1:nrows]:           row index (global)
572                [nrows+1:2*nrows+1]: i-structure index
573     */
574     /*-------------------------------------------*/
575     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
576     buf_si_i    = buf_si + nrows+1;
577     buf_si[0]   = nrows;
578     buf_si_i[0] = 0;
579     nrows       = 0;
580     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
581       nzi = coi[i+1] - coi[i];
582       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
583       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
584       nrows++;
585     }
586     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
587     k++;
588     buf_si += len_si[proc];
589   }
590   for (i=0; i<nrecv; i++) {
591     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
592   }
593   ierr = PetscFree(rwaits);CHKERRQ(ierr);
594   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
595 
596   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
597   ierr = PetscFree(len_ri);CHKERRQ(ierr);
598   ierr = PetscFree(swaits);CHKERRQ(ierr);
599   ierr = PetscFree(buf_s);CHKERRQ(ierr);
600 
601   /* (5) compute the local portion of Cmpi      */
602   /* ------------------------------------------ */
603   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
604   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
605   current_space = free_space;
606 
607   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
608   for (k=0; k<nrecv; k++) {
609     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
610     nrows       = *buf_ri_k[k];
611     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
612     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
613   }
614 
615   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
616   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
617   for (i=0; i<pn; i++) {
618     /* add C_loc into Cmpi */
619     nzi  = c_loc->i[i+1] - c_loc->i[i];
620     Jptr = c_loc->j + c_loc->i[i];
621     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
622 
623     /* add received col data into lnk */
624     for (k=0; k<nrecv; k++) { /* k-th received message */
625       if (i == *nextrow[k]) { /* i-th row */
626         nzi  = *(nextci[k]+1) - *nextci[k];
627         Jptr = buf_rj[k] + *nextci[k];
628         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
629         nextrow[k]++; nextci[k]++;
630       }
631     }
632     nzi = lnk[0];
633 
634     /* copy data into free space, then initialize lnk */
635     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
636     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
637   }
638   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
639   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
640   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
641 
642   /* local sizes and preallocation */
643   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
644   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
645   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
646   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
647 
648   /* members in merge */
649   ierr = PetscFree(id_r);CHKERRQ(ierr);
650   ierr = PetscFree(len_r);CHKERRQ(ierr);
651   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
652   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
653   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
654   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
655   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
656 
657   /* attach the supporting struct to Cmpi for reuse */
658   c = (Mat_MPIAIJ*)Cmpi->data;
659   c->ap           = ptap;
660   ptap->duplicate = Cmpi->ops->duplicate;
661   ptap->destroy   = Cmpi->ops->destroy;
662   ptap->view      = Cmpi->ops->view;
663 
664   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
665   Cmpi->assembled        = PETSC_FALSE;
666   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
667   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
668   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
669   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
670   *C                     = Cmpi;
671 
672    nout = 0;
673    ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_oth->i[ptap->C_oth->rmap->n],c_oth->j,&nout,c_oth->j);CHKERRQ(ierr);
674    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);
675    ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,&nout,c_loc->j);CHKERRQ(ierr);
676    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);
677 
678   PetscFunctionReturn(0);
679 }
680 
681 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
682 {
683   PetscErrorCode      ierr;
684   Mat_APMPI           *ptap;
685   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
686   MPI_Comm            comm;
687   PetscMPIInt         size,rank;
688   Mat                 Cmpi;
689   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
690   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
691   PetscInt            *lnk,i,k,pnz,row,nsend;
692   PetscBT             lnkbt;
693   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
694   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
695   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
696   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
697   MPI_Request         *swaits,*rwaits;
698   MPI_Status          *sstatus,rstatus;
699   PetscLayout         rowmap;
700   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
701   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
702   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
703   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
704   PetscScalar         *apv;
705   PetscTable          ta;
706   MatType             mtype;
707 #if defined(PETSC_USE_INFO)
708   PetscReal           apfill;
709 #endif
710 
711   PetscFunctionBegin;
712   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
713   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
714   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
715 
716   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
717 
718   /* create symbolic parallel matrix Cmpi */
719   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
720   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
721   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
722 
723   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
724   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
725 
726   /* create struct Mat_APMPI and attached it to C later */
727   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
728   ptap->reuse = MAT_INITIAL_MATRIX;
729   ptap->algType = 1;
730 
731   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
732   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
733   /* get P_loc by taking all local rows of P */
734   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
735 
736   /* (0) compute Rd = Pd^T, Ro = Po^T  */
737   /* --------------------------------- */
738   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
739   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
740 
741   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
742   /* ----------------------------------------------------------------- */
743   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
744   if (ptap->P_oth) p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
745 
746   /* create and initialize a linked list */
747   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
748   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
749   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
750   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
751   /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */
752 
753   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
754 
755   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
756   if (ao) {
757     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
758   } else {
759     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
760   }
761   current_space = free_space;
762   nspacedouble  = 0;
763 
764   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
765   api[0] = 0;
766   for (i=0; i<am; i++) {
767     /* diagonal portion: Ad[i,:]*P */
768     ai = ad->i; pi = p_loc->i;
769     nzi = ai[i+1] - ai[i];
770     aj  = ad->j + ai[i];
771     for (j=0; j<nzi; j++) {
772       row  = aj[j];
773       pnz  = pi[row+1] - pi[row];
774       Jptr = p_loc->j + pi[row];
775       /* add non-zero cols of P into the sorted linked list lnk */
776       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
777     }
778     /* off-diagonal portion: Ao[i,:]*P */
779     if (ao) {
780       ai = ao->i; pi = p_oth->i;
781       nzi = ai[i+1] - ai[i];
782       aj  = ao->j + ai[i];
783       for (j=0; j<nzi; j++) {
784         row  = aj[j];
785         pnz  = pi[row+1] - pi[row];
786         Jptr = p_oth->j + pi[row];
787         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
788       }
789     }
790     apnz     = lnk[0];
791     api[i+1] = api[i] + apnz;
792     if (ap_rmax < apnz) ap_rmax = apnz;
793 
794     /* if free space is not available, double the total space in the list */
795     if (current_space->local_remaining<apnz) {
796       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
797       nspacedouble++;
798     }
799 
800     /* Copy data into free space, then initialize lnk */
801     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
802 
803     current_space->array           += apnz;
804     current_space->local_used      += apnz;
805     current_space->local_remaining -= apnz;
806   }
807   /* Allocate space for apj and apv, initialize apj, and */
808   /* destroy list of free space and other temporary array(s) */
809   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
810   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
811   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
812 
813   /* Create AP_loc for reuse */
814   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
815 
816 #if defined(PETSC_USE_INFO)
817   if (ao) {
818     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
819   } else {
820     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
821   }
822   ptap->AP_loc->info.mallocs           = nspacedouble;
823   ptap->AP_loc->info.fill_ratio_given  = fill;
824   ptap->AP_loc->info.fill_ratio_needed = apfill;
825 
826   if (api[am]) {
827     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);
828     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
829   } else {
830     ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
831   }
832 #endif
833 
834   /* (2-1) compute symbolic Co = Ro*AP_loc  */
835   /* ------------------------------------ */
836   ierr = MatSetOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr);
837   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
838 
839   /* (3) send coj of C_oth to other processors  */
840   /* ------------------------------------------ */
841   /* determine row ownership */
842   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
843   rowmap->n  = pn;
844   rowmap->bs = 1;
845   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
846   owners = rowmap->range;
847 
848   /* determine the number of messages to send, their lengths */
849   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
850   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
851   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
852 
853   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
854   coi   = c_oth->i; coj = c_oth->j;
855   con   = ptap->C_oth->rmap->n;
856   proc  = 0;
857   for (i=0; i<con; i++) {
858     while (prmap[i] >= owners[proc+1]) proc++;
859     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
860     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
861   }
862 
863   len          = 0; /* max length of buf_si[], see (4) */
864   owners_co[0] = 0;
865   nsend        = 0;
866   for (proc=0; proc<size; proc++) {
867     owners_co[proc+1] = owners_co[proc] + len_si[proc];
868     if (len_s[proc]) {
869       nsend++;
870       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
871       len         += len_si[proc];
872     }
873   }
874 
875   /* determine the number and length of messages to receive for coi and coj  */
876   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
877   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
878 
879   /* post the Irecv and Isend of coj */
880   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
881   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
882   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
883   for (proc=0, k=0; proc<size; proc++) {
884     if (!len_s[proc]) continue;
885     i    = owners_co[proc];
886     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
887     k++;
888   }
889 
890   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
891   /* ---------------------------------------- */
892   ierr = MatSetOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr);
893   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
894   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
895 
896   /* receives coj are complete */
897   for (i=0; i<nrecv; i++) {
898     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
899   }
900   ierr = PetscFree(rwaits);CHKERRQ(ierr);
901   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
902 
903   /* add received column indices into ta to update Crmax */
904   for (k=0; k<nrecv; k++) {/* k-th received message */
905     Jptr = buf_rj[k];
906     for (j=0; j<len_r[k]; j++) {
907       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
908     }
909   }
910   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
911   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
912 
913   /* (4) send and recv coi */
914   /*-----------------------*/
915   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
916   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
917   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
918   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
919   for (proc=0,k=0; proc<size; proc++) {
920     if (!len_s[proc]) continue;
921     /* form outgoing message for i-structure:
922          buf_si[0]:                 nrows to be sent
923                [1:nrows]:           row index (global)
924                [nrows+1:2*nrows+1]: i-structure index
925     */
926     /*-------------------------------------------*/
927     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
928     buf_si_i    = buf_si + nrows+1;
929     buf_si[0]   = nrows;
930     buf_si_i[0] = 0;
931     nrows       = 0;
932     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
933       nzi = coi[i+1] - coi[i];
934       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
935       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
936       nrows++;
937     }
938     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
939     k++;
940     buf_si += len_si[proc];
941   }
942   for (i=0; i<nrecv; i++) {
943     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
944   }
945   ierr = PetscFree(rwaits);CHKERRQ(ierr);
946   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
947 
948   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
949   ierr = PetscFree(len_ri);CHKERRQ(ierr);
950   ierr = PetscFree(swaits);CHKERRQ(ierr);
951   ierr = PetscFree(buf_s);CHKERRQ(ierr);
952 
953   /* (5) compute the local portion of Cmpi      */
954   /* ------------------------------------------ */
955   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
956   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
957   current_space = free_space;
958 
959   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
960   for (k=0; k<nrecv; k++) {
961     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
962     nrows       = *buf_ri_k[k];
963     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
964     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
965   }
966 
967   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
968   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
969   for (i=0; i<pn; i++) {
970     /* add C_loc into Cmpi */
971     nzi  = c_loc->i[i+1] - c_loc->i[i];
972     Jptr = c_loc->j + c_loc->i[i];
973     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
974 
975     /* add received col data into lnk */
976     for (k=0; k<nrecv; k++) { /* k-th received message */
977       if (i == *nextrow[k]) { /* i-th row */
978         nzi  = *(nextci[k]+1) - *nextci[k];
979         Jptr = buf_rj[k] + *nextci[k];
980         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
981         nextrow[k]++; nextci[k]++;
982       }
983     }
984     nzi = lnk[0];
985 
986     /* copy data into free space, then initialize lnk */
987     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
988     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
989   }
990   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
991   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
992   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
993 
994   /* local sizes and preallocation */
995   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
996   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
997   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
998   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
999 
1000   /* members in merge */
1001   ierr = PetscFree(id_r);CHKERRQ(ierr);
1002   ierr = PetscFree(len_r);CHKERRQ(ierr);
1003   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
1004   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
1005   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
1006   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
1007   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
1008 
1009   /* attach the supporting struct to Cmpi for reuse */
1010   c = (Mat_MPIAIJ*)Cmpi->data;
1011   c->ap           = ptap;
1012   ptap->duplicate = Cmpi->ops->duplicate;
1013   ptap->destroy   = Cmpi->ops->destroy;
1014   ptap->view      = Cmpi->ops->view;
1015   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
1016 
1017   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1018   Cmpi->assembled        = PETSC_FALSE;
1019   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1020   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
1021   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1022   *C                     = Cmpi;
1023   PetscFunctionReturn(0);
1024 }
1025 
1026 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
1027 {
1028   PetscErrorCode    ierr;
1029   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1030   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1031   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
1032   Mat_APMPI         *ptap = c->ap;
1033   Mat               AP_loc,C_loc,C_oth;
1034   PetscInt          i,rstart,rend,cm,ncols,row;
1035   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
1036   PetscScalar       *apa;
1037   const PetscInt    *cols;
1038   const PetscScalar *vals;
1039 
1040   PetscFunctionBegin;
1041   if (!ptap) {
1042     MPI_Comm comm;
1043     ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1044     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1045   }
1046 
1047   ierr = MatZeroEntries(C);CHKERRQ(ierr);
1048   /* 1) get R = Pd^T,Ro = Po^T */
1049   if (ptap->reuse == MAT_REUSE_MATRIX) {
1050     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1051     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1052   }
1053 
1054   /* 2) get AP_loc */
1055   AP_loc = ptap->AP_loc;
1056   ap = (Mat_SeqAIJ*)AP_loc->data;
1057 
1058   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
1059   /*-----------------------------------------------------*/
1060   if (ptap->reuse == MAT_REUSE_MATRIX) {
1061     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
1062     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1063     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1064   }
1065 
1066   /* 2-2) compute numeric A_loc*P - dominating part */
1067   /* ---------------------------------------------- */
1068   /* get data from symbolic products */
1069   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1070   if (ptap->P_oth) {
1071     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1072   }
1073   apa   = ptap->apa;
1074   api   = ap->i;
1075   apj   = ap->j;
1076   for (i=0; i<am; i++) {
1077     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
1078     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
1079     apnz = api[i+1] - api[i];
1080     for (j=0; j<apnz; j++) {
1081       col = apj[j+api[i]];
1082       ap->a[j+ap->i[i]] = apa[col];
1083       apa[col] = 0.0;
1084     }
1085     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1086   }
1087 
1088   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
1089   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
1090   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
1091   C_loc = ptap->C_loc;
1092   C_oth = ptap->C_oth;
1093 
1094   /* add C_loc and Co to to C */
1095   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
1096 
1097   /* C_loc -> C */
1098   cm    = C_loc->rmap->N;
1099   c_seq = (Mat_SeqAIJ*)C_loc->data;
1100   cols = c_seq->j;
1101   vals = c_seq->a;
1102 
1103 
1104   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
1105   /* when there are no off-processor parts.  */
1106   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
1107   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
1108   /* a table, and other, more complex stuff has to be done. */
1109   if (C->assembled) {
1110     C->was_assembled = PETSC_TRUE;
1111     C->assembled     = PETSC_FALSE;
1112   }
1113   if (C->was_assembled) {
1114     for (i=0; i<cm; i++) {
1115       ncols = c_seq->i[i+1] - c_seq->i[i];
1116       row = rstart + i;
1117       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1118       cols += ncols; vals += ncols;
1119     }
1120   } else {
1121     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
1122   }
1123 
1124   /* Co -> C, off-processor part */
1125   cm = C_oth->rmap->N;
1126   c_seq = (Mat_SeqAIJ*)C_oth->data;
1127   cols = c_seq->j;
1128   vals = c_seq->a;
1129   for (i=0; i<cm; i++) {
1130     ncols = c_seq->i[i+1] - c_seq->i[i];
1131     row = p->garray[i];
1132     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1133     cols += ncols; vals += ncols;
1134   }
1135 
1136   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1137   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1138 
1139   ptap->reuse = MAT_REUSE_MATRIX;
1140 
1141   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
1142   if (ptap->freestruct) {
1143     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
1144   }
1145   PetscFunctionReturn(0);
1146 }
1147