xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 580bdb303e1ee3b1222b2042810b4c26340259c6)
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 #include <petsc/private/hashmapiv.h>
13 #include <petsc/private/hashseti.h>
14 #include <petscsf.h>
15 
16 
17 PetscErrorCode MatView_MPIAIJ_PtAP(Mat A,PetscViewer viewer)
18 {
19   PetscErrorCode    ierr;
20   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data;
21   Mat_APMPI         *ptap=a->ap;
22   PetscBool         iascii;
23   PetscViewerFormat format;
24 
25   PetscFunctionBegin;
26   if (!ptap) {
27     /* hack: MatDuplicate() sets oldmat->ops->view to newmat which is a base mat class with null ptpa! */
28     A->ops->view = MatView_MPIAIJ;
29     ierr = (A->ops->view)(A,viewer);CHKERRQ(ierr);
30     PetscFunctionReturn(0);
31   }
32 
33   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
34   if (iascii) {
35     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
36     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
37       if (ptap->algType == 0) {
38         ierr = PetscViewerASCIIPrintf(viewer,"using scalable MatPtAP() implementation\n");CHKERRQ(ierr);
39       } else if (ptap->algType == 1) {
40         ierr = PetscViewerASCIIPrintf(viewer,"using nonscalable MatPtAP() implementation\n");CHKERRQ(ierr);
41       } else if (ptap->algType == 2) {
42         ierr = PetscViewerASCIIPrintf(viewer,"using allatonce MatPtAP() implementation\n");CHKERRQ(ierr);
43       } else if (ptap->algType == 3) {
44         ierr = PetscViewerASCIIPrintf(viewer,"using merged allatonce MatPtAP() implementation\n");CHKERRQ(ierr);
45       }
46     }
47   }
48   ierr = (ptap->view)(A,viewer);CHKERRQ(ierr);
49   PetscFunctionReturn(0);
50 }
51 
52 PetscErrorCode MatFreeIntermediateDataStructures_MPIAIJ_AP(Mat A)
53 {
54   PetscErrorCode      ierr;
55   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data;
56   Mat_APMPI           *ptap=a->ap;
57   Mat_Merge_SeqsToMPI *merge;
58 
59   PetscFunctionBegin;
60   if (!ptap) PetscFunctionReturn(0);
61 
62   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
63   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
64   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
65   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
66   ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
67   ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
68   ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
69   if (ptap->AP_loc) { /* used by alg_rap */
70     Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
71     ierr = PetscFree(ap->i);CHKERRQ(ierr);
72     ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
73     ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
74   } else { /* used by alg_ptap */
75     ierr = PetscFree(ptap->api);CHKERRQ(ierr);
76     ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
77   }
78   ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
79   ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
80   if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
81 
82   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
83 
84   merge=ptap->merge;
85   if (merge) { /* used by alg_ptap */
86     ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
87     ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
88     ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
89     ierr = PetscFree(merge->bi);CHKERRQ(ierr);
90     ierr = PetscFree(merge->bj);CHKERRQ(ierr);
91     ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
92     ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
93     ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
94     ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
95     ierr = PetscFree(merge->coi);CHKERRQ(ierr);
96     ierr = PetscFree(merge->coj);CHKERRQ(ierr);
97     ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
98     ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
99     ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
100   }
101   ierr = ISLocalToGlobalMappingDestroy(&ptap->ltog);CHKERRQ(ierr);
102 
103   ierr = PetscSFDestroy(&ptap->sf);CHKERRQ(ierr);
104   ierr = PetscFree(ptap->c_othi);CHKERRQ(ierr);
105   ierr = PetscFree(ptap->c_rmti);CHKERRQ(ierr);
106   PetscFunctionReturn(0);
107 }
108 
109 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
110 {
111   PetscErrorCode ierr;
112   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
113   Mat_APMPI      *ptap=a->ap;
114 
115   PetscFunctionBegin;
116   ierr = (*A->ops->freeintermediatedatastructures)(A);CHKERRQ(ierr);
117   ierr = (*ptap->destroy)(A);CHKERRQ(ierr); /* MatDestroy_MPIAIJ(A) */
118   ierr = PetscFree(ptap);CHKERRQ(ierr);
119   PetscFunctionReturn(0);
120 }
121 
122 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
123 {
124   PetscErrorCode ierr;
125   PetscBool      flg;
126   MPI_Comm       comm;
127 #if !defined(PETSC_HAVE_HYPRE)
128   const char          *algTypes[4] = {"scalable","nonscalable","allatonce","allatonce_merged"};
129   PetscInt            nalg=4;
130 #else
131   const char          *algTypes[5] = {"scalable","nonscalable","allatonce","allatonce_merged","hypre"};
132   PetscInt            nalg=5;
133 #endif
134   PetscInt            pN=P->cmap->N,alg=1; /* set default algorithm */
135 
136   PetscFunctionBegin;
137   /* check if matrix local sizes are compatible */
138   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
139   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);
140   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);
141 
142   if (scall == MAT_INITIAL_MATRIX) {
143     /* pick an algorithm */
144     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
145     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
146     ierr = PetscOptionsEnd();CHKERRQ(ierr);
147 
148     if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */
149       MatInfo     Ainfo,Pinfo;
150       PetscInt    nz_local;
151       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
152 
153       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
154       ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
155       nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
156 
157       if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
158       ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
159 
160       if (alg_scalable) {
161         alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
162       }
163     }
164 
165     switch (alg) {
166     case 1:
167       /* do R=P^T locally, then C=R*A*P -- nonscalable */
168       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
169       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
170       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
171       break;
172     case 2:
173       /* compute C=P^T*A*P allatonce */
174       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
175       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(A,P,fill,C);CHKERRQ(ierr);
176       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
177       break;
178     case 3:
179       /* compute C=P^T*A*P allatonce */
180       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
181       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(A,P,fill,C);CHKERRQ(ierr);
182       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
183       break;
184 #if defined(PETSC_HAVE_HYPRE)
185     case 4:
186       /* Use boomerAMGBuildCoarseOperator */
187       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
188       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
189       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
190       break;
191 #endif
192     default:
193       /* do R=P^T locally, then C=R*A*P */
194       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
195       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr);
196       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
197       break;
198     }
199 
200     if (alg == 0 || alg == 1 || alg == 2 || alg == 3) {
201       Mat_MPIAIJ *c  = (Mat_MPIAIJ*)(*C)->data;
202       Mat_APMPI  *ap = c->ap;
203       ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr);
204       ap->freestruct = PETSC_FALSE;
205       ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr);
206       ierr = PetscOptionsEnd();CHKERRQ(ierr);
207     }
208   }
209 
210   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
211   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
212   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
213   PetscFunctionReturn(0);
214 }
215 
216 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C)
217 {
218   PetscErrorCode    ierr;
219   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
220   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
221   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
222   Mat_APMPI         *ptap = c->ap;
223   Mat               AP_loc,C_loc,C_oth;
224   PetscInt          i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz,nout;
225   PetscScalar       *apa;
226   const PetscInt    *cols;
227   const PetscScalar *vals;
228 
229   PetscFunctionBegin;
230   if (!ptap->AP_loc) {
231     MPI_Comm comm;
232     ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
233     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
234   }
235 
236   ierr = MatZeroEntries(C);CHKERRQ(ierr);
237 
238   /* 1) get R = Pd^T,Ro = Po^T */
239   if (ptap->reuse == MAT_REUSE_MATRIX) {
240     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
241     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
242   }
243 
244   /* 2) get AP_loc */
245   AP_loc = ptap->AP_loc;
246   ap = (Mat_SeqAIJ*)AP_loc->data;
247 
248   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
249   /*-----------------------------------------------------*/
250   if (ptap->reuse == MAT_REUSE_MATRIX) {
251     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
252     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
253     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
254   }
255 
256   /* 2-2) compute numeric A_loc*P - dominating part */
257   /* ---------------------------------------------- */
258   /* get data from symbolic products */
259   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
260   if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
261 
262   api   = ap->i;
263   apj   = ap->j;
264   ierr = ISLocalToGlobalMappingApply(ptap->ltog,api[AP_loc->rmap->n],apj,apj);CHKERRQ(ierr);
265   for (i=0; i<am; i++) {
266     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
267     apnz = api[i+1] - api[i];
268     apa = ap->a + api[i];
269     ierr = PetscArrayzero(apa,apnz);CHKERRQ(ierr);
270     AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa);
271   }
272   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,api[AP_loc->rmap->n],apj,&nout,apj);CHKERRQ(ierr);
273   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);
274 
275   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
276   /* Always use scalable version since we are in the MPI scalable version */
277   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
278   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
279 
280   C_loc = ptap->C_loc;
281   C_oth = ptap->C_oth;
282 
283   /* add C_loc and Co to to C */
284   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
285 
286   /* C_loc -> C */
287   cm    = C_loc->rmap->N;
288   c_seq = (Mat_SeqAIJ*)C_loc->data;
289   cols = c_seq->j;
290   vals = c_seq->a;
291   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_loc->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr);
292 
293   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
294   /* when there are no off-processor parts.  */
295   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
296   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
297   /* a table, and other, more complex stuff has to be done. */
298   if (C->assembled) {
299     C->was_assembled = PETSC_TRUE;
300     C->assembled     = PETSC_FALSE;
301   }
302   if (C->was_assembled) {
303     for (i=0; i<cm; i++) {
304       ncols = c_seq->i[i+1] - c_seq->i[i];
305       row = rstart + i;
306       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
307       cols += ncols; vals += ncols;
308     }
309   } else {
310     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
311   }
312   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_loc->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr);
313   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);
314 
315   /* Co -> C, off-processor part */
316   cm = C_oth->rmap->N;
317   c_seq = (Mat_SeqAIJ*)C_oth->data;
318   cols = c_seq->j;
319   vals = c_seq->a;
320   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_oth->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr);
321   for (i=0; i<cm; i++) {
322     ncols = c_seq->i[i+1] - c_seq->i[i];
323     row = p->garray[i];
324     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
325     cols += ncols; vals += ncols;
326   }
327   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
328   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
329 
330   ptap->reuse = MAT_REUSE_MATRIX;
331 
332   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_oth->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr);
333   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);
334 
335   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
336   if (ptap->freestruct) {
337     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
338   }
339   PetscFunctionReturn(0);
340 }
341 
342 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
343 {
344   PetscErrorCode      ierr;
345   Mat_APMPI           *ptap;
346   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
347   MPI_Comm            comm;
348   PetscMPIInt         size,rank;
349   Mat                 Cmpi,P_loc,P_oth;
350   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
351   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
352   PetscInt            *lnk,i,k,pnz,row,nsend;
353   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
354   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
355   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
356   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
357   MPI_Request         *swaits,*rwaits;
358   MPI_Status          *sstatus,rstatus;
359   PetscLayout         rowmap;
360   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
361   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
362   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi,nout;
363   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
364   PetscScalar         *apv;
365   PetscTable          ta;
366   MatType             mtype;
367   const char          *prefix;
368 #if defined(PETSC_USE_INFO)
369   PetscReal           apfill;
370 #endif
371 
372   PetscFunctionBegin;
373   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
374   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
375   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
376 
377   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
378 
379   /* create symbolic parallel matrix Cmpi */
380   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
381   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
382   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
383 
384   /* create struct Mat_APMPI and attached it to C later */
385   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
386   ptap->reuse = MAT_INITIAL_MATRIX;
387   ptap->algType = 0;
388 
389   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
390   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
391   /* get P_loc by taking all local rows of P */
392   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
393 
394   ptap->P_loc = P_loc;
395   ptap->P_oth = P_oth;
396 
397   /* (0) compute Rd = Pd^T, Ro = Po^T  */
398   /* --------------------------------- */
399   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
400   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
401 
402   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
403   /* ----------------------------------------------------------------- */
404   p_loc  = (Mat_SeqAIJ*)P_loc->data;
405   if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data;
406 
407   /* create and initialize a linked list */
408   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
409   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
410   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
411   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
412 
413   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
414 
415   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
416   if (ao) {
417     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
418   } else {
419     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
420   }
421   current_space = free_space;
422   nspacedouble  = 0;
423 
424   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
425   api[0] = 0;
426   for (i=0; i<am; i++) {
427     /* diagonal portion: Ad[i,:]*P */
428     ai = ad->i; pi = p_loc->i;
429     nzi = ai[i+1] - ai[i];
430     aj  = ad->j + ai[i];
431     for (j=0; j<nzi; j++) {
432       row  = aj[j];
433       pnz  = pi[row+1] - pi[row];
434       Jptr = p_loc->j + pi[row];
435       /* add non-zero cols of P into the sorted linked list lnk */
436       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
437     }
438     /* off-diagonal portion: Ao[i,:]*P */
439     if (ao) {
440       ai = ao->i; pi = p_oth->i;
441       nzi = ai[i+1] - ai[i];
442       aj  = ao->j + ai[i];
443       for (j=0; j<nzi; j++) {
444         row  = aj[j];
445         pnz  = pi[row+1] - pi[row];
446         Jptr = p_oth->j + pi[row];
447         ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
448       }
449     }
450     apnz     = lnk[0];
451     api[i+1] = api[i] + apnz;
452 
453     /* if free space is not available, double the total space in the list */
454     if (current_space->local_remaining<apnz) {
455       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
456       nspacedouble++;
457     }
458 
459     /* Copy data into free space, then initialize lnk */
460     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
461 
462     current_space->array           += apnz;
463     current_space->local_used      += apnz;
464     current_space->local_remaining -= apnz;
465   }
466   /* Allocate space for apj and apv, initialize apj, and */
467   /* destroy list of free space and other temporary array(s) */
468   ierr = PetscCalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
469   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
470   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
471 
472   /* Create AP_loc for reuse */
473   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
474   ierr = MatSeqAIJCompactOutExtraColumns_SeqAIJ(ptap->AP_loc, &ptap->ltog);CHKERRQ(ierr);
475 
476 #if defined(PETSC_USE_INFO)
477   if (ao) {
478     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
479   } else {
480     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
481   }
482   ptap->AP_loc->info.mallocs           = nspacedouble;
483   ptap->AP_loc->info.fill_ratio_given  = fill;
484   ptap->AP_loc->info.fill_ratio_needed = apfill;
485 
486   if (api[am]) {
487     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);
488     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
489   } else {
490     ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
491   }
492 #endif
493 
494   /* (2-1) compute symbolic Co = Ro*AP_loc  */
495   /* ------------------------------------ */
496   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
497   ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr);
498   ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr);
499   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
500 
501   /* (3) send coj of C_oth to other processors  */
502   /* ------------------------------------------ */
503   /* determine row ownership */
504   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
505   rowmap->n  = pn;
506   rowmap->bs = 1;
507   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
508   owners = rowmap->range;
509 
510   /* determine the number of messages to send, their lengths */
511   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
512   /* TODO: Use PetscCalloc2() */
513   ierr = PetscArrayzero(len_s,size);CHKERRQ(ierr);
514   ierr = PetscArrayzero(len_si,size);CHKERRQ(ierr);
515 
516   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
517   coi   = c_oth->i; coj = c_oth->j;
518   con   = ptap->C_oth->rmap->n;
519   proc  = 0;
520   ierr = ISLocalToGlobalMappingApply(ptap->ltog,coi[con],coj,coj);CHKERRQ(ierr);
521   for (i=0; i<con; i++) {
522     while (prmap[i] >= owners[proc+1]) proc++;
523     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
524     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
525   }
526 
527   len          = 0; /* max length of buf_si[], see (4) */
528   owners_co[0] = 0;
529   nsend        = 0;
530   for (proc=0; proc<size; proc++) {
531     owners_co[proc+1] = owners_co[proc] + len_si[proc];
532     if (len_s[proc]) {
533       nsend++;
534       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
535       len         += len_si[proc];
536     }
537   }
538 
539   /* determine the number and length of messages to receive for coi and coj  */
540   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
541   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
542 
543   /* post the Irecv and Isend of coj */
544   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
545   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
546   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
547   for (proc=0, k=0; proc<size; proc++) {
548     if (!len_s[proc]) continue;
549     i    = owners_co[proc];
550     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
551     k++;
552   }
553 
554   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
555   /* ---------------------------------------- */
556   ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr);
557   ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr);
558   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
559   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
560   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,c_loc->j);CHKERRQ(ierr);
561 
562   /* receives coj are complete */
563   for (i=0; i<nrecv; i++) {
564     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
565   }
566   ierr = PetscFree(rwaits);CHKERRQ(ierr);
567   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
568 
569   /* add received column indices into ta to update Crmax */
570   for (k=0; k<nrecv; k++) {/* k-th received message */
571     Jptr = buf_rj[k];
572     for (j=0; j<len_r[k]; j++) {
573       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
574     }
575   }
576   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
577   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
578 
579   /* (4) send and recv coi */
580   /*-----------------------*/
581   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
582   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
583   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
584   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
585   for (proc=0,k=0; proc<size; proc++) {
586     if (!len_s[proc]) continue;
587     /* form outgoing message for i-structure:
588          buf_si[0]:                 nrows to be sent
589                [1:nrows]:           row index (global)
590                [nrows+1:2*nrows+1]: i-structure index
591     */
592     /*-------------------------------------------*/
593     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
594     buf_si_i    = buf_si + nrows+1;
595     buf_si[0]   = nrows;
596     buf_si_i[0] = 0;
597     nrows       = 0;
598     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
599       nzi = coi[i+1] - coi[i];
600       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
601       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
602       nrows++;
603     }
604     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
605     k++;
606     buf_si += len_si[proc];
607   }
608   for (i=0; i<nrecv; i++) {
609     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
610   }
611   ierr = PetscFree(rwaits);CHKERRQ(ierr);
612   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
613 
614   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
615   ierr = PetscFree(len_ri);CHKERRQ(ierr);
616   ierr = PetscFree(swaits);CHKERRQ(ierr);
617   ierr = PetscFree(buf_s);CHKERRQ(ierr);
618 
619   /* (5) compute the local portion of Cmpi      */
620   /* ------------------------------------------ */
621   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
622   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
623   current_space = free_space;
624 
625   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
626   for (k=0; k<nrecv; k++) {
627     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
628     nrows       = *buf_ri_k[k];
629     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
630     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
631   }
632 
633   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
634   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
635   for (i=0; i<pn; i++) {
636     /* add C_loc into Cmpi */
637     nzi  = c_loc->i[i+1] - c_loc->i[i];
638     Jptr = c_loc->j + c_loc->i[i];
639     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
640 
641     /* add received col data into lnk */
642     for (k=0; k<nrecv; k++) { /* k-th received message */
643       if (i == *nextrow[k]) { /* i-th row */
644         nzi  = *(nextci[k]+1) - *nextci[k];
645         Jptr = buf_rj[k] + *nextci[k];
646         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
647         nextrow[k]++; nextci[k]++;
648       }
649     }
650     nzi = lnk[0];
651 
652     /* copy data into free space, then initialize lnk */
653     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
654     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
655   }
656   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
657   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
658   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
659 
660   /* local sizes and preallocation */
661   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
662   if (P->cmap->bs > 0) {
663     ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr);
664     ierr = PetscLayoutSetBlockSize(Cmpi->cmap,P->cmap->bs);CHKERRQ(ierr);
665   }
666   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
667   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
668 
669   /* members in merge */
670   ierr = PetscFree(id_r);CHKERRQ(ierr);
671   ierr = PetscFree(len_r);CHKERRQ(ierr);
672   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
673   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
674   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
675   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
676   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
677 
678   /* attach the supporting struct to Cmpi for reuse */
679   c = (Mat_MPIAIJ*)Cmpi->data;
680   c->ap           = ptap;
681   ptap->duplicate = Cmpi->ops->duplicate;
682   ptap->destroy   = Cmpi->ops->destroy;
683   ptap->view      = Cmpi->ops->view;
684 
685   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
686   Cmpi->assembled        = PETSC_FALSE;
687   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
688   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
689   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
690   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
691   *C                     = Cmpi;
692 
693    nout = 0;
694    ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_oth->i[ptap->C_oth->rmap->n],c_oth->j,&nout,c_oth->j);CHKERRQ(ierr);
695    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);
696    ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,&nout,c_loc->j);CHKERRQ(ierr);
697    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);
698 
699   PetscFunctionReturn(0);
700 }
701 
702 PETSC_STATIC_INLINE PetscErrorCode MatPtAPSymbolicComputeOneRowOfAP_private(Mat A,Mat P,Mat P_oth,PetscInt i,PetscHSetI dht,PetscHSetI oht)
703 {
704   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
705   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_oth=(Mat_SeqAIJ*)P_oth->data,*pd=(Mat_SeqAIJ*)p->A->data,*po=(Mat_SeqAIJ*)p->B->data;
706   PetscInt             *ai,nzi,j,*aj,row,col,*pi,*pj,pnz,nzpi,*p_othcols,k;
707   PetscInt             pcstart,pcend,column;
708   PetscErrorCode       ierr;
709 
710   PetscFunctionBegin;
711 
712   pcstart = P->cmap->rstart;
713   pcend   = P->cmap->rend;
714   /* diagonal portion: Ad[i,:]*P */
715   ai = ad->i;
716   nzi = ai[i+1] - ai[i];
717   aj  = ad->j + ai[i];
718   for (j=0; j<nzi; j++) {
719     row  = aj[j];
720     nzpi = pd->i[row+1] - pd->i[row];
721     pj  = pd->j + pd->i[row];
722     for (k=0; k<nzpi; k++) {
723       ierr = PetscHSetIAdd(dht,pj[k]+pcstart);CHKERRQ(ierr);
724     }
725   }
726   for (j=0; j<nzi; j++) {
727     row  = aj[j];
728     nzpi = po->i[row+1] - po->i[row];
729     pj  = po->j + po->i[row];
730     for (k=0; k<nzpi; k++) {
731       ierr = PetscHSetIAdd(oht,p->garray[pj[k]]);CHKERRQ(ierr);
732     }
733   }
734 
735   /* off diagonal part: Ao[i, :]*P_oth */
736   if (ao) {
737     ai = ao->i;
738     pi = p_oth->i;
739     nzi = ai[i+1] - ai[i];
740     aj  = ao->j + ai[i];
741     for (j=0; j<nzi; j++) {
742       row  = aj[j];
743       pnz  = pi[row+1] - pi[row];
744       p_othcols = p_oth->j + pi[row];
745       for (col=0; col<pnz; col++) {
746         column = p_othcols[col];
747         if (column>=pcstart && column<pcend) {
748           ierr = PetscHSetIAdd(dht,column);CHKERRQ(ierr);
749         } else {
750           ierr = PetscHSetIAdd(oht,column);CHKERRQ(ierr);
751         }
752       }
753     }
754   } /* end if (ao) */
755   PetscFunctionReturn(0);
756 }
757 
758 PETSC_STATIC_INLINE PetscErrorCode MatPtAPNumericComputeOneRowOfAP_private(Mat A,Mat P,Mat P_oth,PetscInt i,PetscHMapIV hmap)
759 {
760   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
761   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_oth=(Mat_SeqAIJ*)P_oth->data,*pd=(Mat_SeqAIJ*)p->A->data,*po=(Mat_SeqAIJ*)p->B->data;
762   PetscInt             *ai,nzi,j,*aj,row,col,*pi,pnz,*p_othcols,pcstart,*pj,k,nzpi;
763   PetscScalar          ra,*aa,*pa;
764   PetscErrorCode       ierr;
765 
766   PetscFunctionBegin;
767   pcstart = P->cmap->rstart;
768   /* diagonal portion: Ad[i,:]*P */
769   ai  = ad->i;
770   nzi = ai[i+1] - ai[i];
771   aj  = ad->j + ai[i];
772   aa  = ad->a + ai[i];
773 
774   for (j=0; j<nzi; j++) {
775     ra   = aa[j];
776     row  = aj[j];
777     nzpi = pd->i[row+1] - pd->i[row];
778     pj = pd->j + pd->i[row];
779     pa = pd->a + pd->i[row];
780     for (k=0; k<nzpi; k++) {
781       ierr = PetscHMapIVAddValue(hmap,pj[k]+pcstart,ra*pa[k]);CHKERRQ(ierr);
782     }
783     ierr = PetscLogFlops(2.0*nzpi);CHKERRQ(ierr);
784   }
785   for (j=0; j<nzi; j++) {
786     ra   = aa[j];
787     row  = aj[j];
788     nzpi = po->i[row+1] - po->i[row];
789     pj = po->j + po->i[row];
790     pa = po->a + po->i[row];
791     for (k=0; k<nzpi; k++) {
792       ierr = PetscHMapIVAddValue(hmap,p->garray[pj[k]],ra*pa[k]);CHKERRQ(ierr);
793     }
794     ierr = PetscLogFlops(2.0*nzpi);CHKERRQ(ierr);
795   }
796 
797 
798   /* off diagonal part: Ao[i, :]*P_oth */
799   if (ao) {
800     ai = ao->i;
801     pi = p_oth->i;
802     nzi = ai[i+1] - ai[i];
803     aj  = ao->j + ai[i];
804     aa  = ao->a + ai[i];
805     for (j=0; j<nzi; j++) {
806       row  = aj[j];
807       ra   = aa[j];
808       pnz  = pi[row+1] - pi[row];
809       p_othcols = p_oth->j + pi[row];
810       pa   = p_oth->a + pi[row];
811       for (col=0; col<pnz; col++) {
812         ierr = PetscHMapIVAddValue(hmap,p_othcols[col],ra*pa[col]);CHKERRQ(ierr);
813       }
814       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
815     }
816   } /* end if (ao) */
817 
818   PetscFunctionReturn(0);
819 }
820 
821 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,Mat C)
822 {
823   PetscErrorCode    ierr;
824   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
825   Mat_SeqAIJ        *cd,*co,*po,*pd;
826   Mat_APMPI         *ptap = c->ap;
827   PetscHMapIV       hmap;
828   PetscInt          i,j,jj,kk,nzi,*c_rmtj,voff,*c_othj,pn,pon,pcstart,pcend,ccstart,ccend,row,am,*poj,*pdj,*apindices,cmaxr,*c_rmtc,*c_rmtjj,*dcc,*occ,loc;
829   PetscScalar       *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa;
830   MPI_Comm          comm;
831 
832   PetscFunctionBegin;
833   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
834   if (!ptap->P_oth) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
835 
836   ierr = MatZeroEntries(C);CHKERRQ(ierr);
837 
838   /* Get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
839   /*-----------------------------------------------------*/
840   if (ptap->reuse == MAT_REUSE_MATRIX) {
841     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
842     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
843   }
844 
845   po = (Mat_SeqAIJ*) p->B->data;
846   pd = (Mat_SeqAIJ*) p->A->data;
847 
848   ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr);
849 
850   ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr);
851   ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr);
852   cmaxr = 0;
853   for (i=0; i<pon; i++) {
854     cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]);
855   }
856   ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr);
857   ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr);
858   ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr);
859   for (i=0; i<am && pon; i++) {
860     ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr);
861     nzi = po->i[i+1] - po->i[i];
862     if (!nzi) continue;
863     ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr);
864     voff = 0;
865     ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr);
866     if (!voff) continue;
867 
868     /* Form C(ii, :) */
869     poj = po->j + po->i[i];
870     poa = po->a + po->i[i];
871     for (j=0; j<nzi; j++) {
872       c_rmtjj = c_rmtj + ptap->c_rmti[poj[j]];
873       c_rmtaa = c_rmta + ptap->c_rmti[poj[j]];
874       for (jj=0; jj<voff; jj++) {
875         apvaluestmp[jj] = apvalues[jj]*poa[j];
876         /*If the row is empty */
877         if (!c_rmtc[poj[j]]) {
878           c_rmtjj[jj] = apindices[jj];
879           c_rmtaa[jj] = apvaluestmp[jj];
880           c_rmtc[poj[j]]++;
881         } else {
882           ierr = PetscFindInt(apindices[jj],c_rmtc[poj[j]],c_rmtjj,&loc);CHKERRQ(ierr);
883           if (loc>=0){ /* hit */
884             c_rmtaa[loc] += apvaluestmp[jj];
885             ierr = PetscLogFlops(1.0);CHKERRQ(ierr);
886           } else { /* new element */
887             loc = -(loc+1);
888             /* Move data backward */
889             for (kk=c_rmtc[poj[j]]; kk>loc; kk--) {
890               c_rmtjj[kk] = c_rmtjj[kk-1];
891               c_rmtaa[kk] = c_rmtaa[kk-1];
892             }/* End kk */
893             c_rmtjj[loc] = apindices[jj];
894             c_rmtaa[loc] = apvaluestmp[jj];
895             c_rmtc[poj[j]]++;
896           }
897         }
898         ierr = PetscLogFlops(voff);CHKERRQ(ierr);
899       } /* End jj */
900     } /* End j */
901   } /* End i */
902 
903   ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr);
904 
905   ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr);
906   ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr);
907 
908   ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
909   ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr);
910   ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr);
911   cd = (Mat_SeqAIJ*)(c->A)->data;
912   co = (Mat_SeqAIJ*)(c->B)->data;
913 
914   cmaxr = 0;
915   for (i=0; i<pn; i++) {
916     cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i]));
917   }
918   ierr = PetscCalloc5(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pn,&dcc,pn,&occ);CHKERRQ(ierr);
919   ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr);
920   ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr);
921   for (i=0; i<am && pn; i++) {
922     ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr);
923     nzi = pd->i[i+1] - pd->i[i];
924     if (!nzi) continue;
925     ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr);
926     voff = 0;
927     ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr);
928     if (!voff) continue;
929     /* Form C(ii, :) */
930     pdj = pd->j + pd->i[i];
931     pda = pd->a + pd->i[i];
932     for (j=0; j<nzi; j++) {
933       row = pcstart + pdj[j];
934       for (jj=0; jj<voff; jj++) {
935         apvaluestmp[jj] = apvalues[jj]*pda[j];
936       }
937       ierr = PetscLogFlops(voff);CHKERRQ(ierr);
938       ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr);
939     }
940   }
941 
942   ierr = MatGetOwnershipRangeColumn(C,&ccstart,&ccend);CHKERRQ(ierr);
943   ierr = PetscFree5(apindices,apvalues,apvaluestmp,dcc,occ);CHKERRQ(ierr);
944   ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr);
945   ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_SUM);CHKERRQ(ierr);
946   ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr);
947   ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr);
948 
949   /* Add contributions from remote */
950   for (i = 0; i < pn; i++) {
951     row = i + pcstart;
952     ierr = MatSetValues(C,1,&row,ptap->c_othi[i+1]-ptap->c_othi[i],c_othj+ptap->c_othi[i],c_otha+ptap->c_othi[i],ADD_VALUES);CHKERRQ(ierr);
953   }
954   ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr);
955 
956   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
957   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
958 
959   ptap->reuse = MAT_REUSE_MATRIX;
960 
961   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
962   if (ptap->freestruct) {
963     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
964   }
965   PetscFunctionReturn(0);
966 }
967 
968 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,Mat C)
969 {
970   PetscErrorCode    ierr;
971   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
972   Mat_SeqAIJ        *cd,*co,*po,*pd;
973   Mat_APMPI         *ptap = c->ap;
974   PetscHMapIV       hmap;
975   PetscInt          i,j,jj,kk,nzi,dnzi,*c_rmtj,voff,*c_othj,pn,pon,pcstart,pcend,row,am,*poj,*pdj,*apindices,cmaxr,*c_rmtc,*c_rmtjj,loc;
976   PetscScalar       *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa;
977   MPI_Comm          comm;
978 
979   PetscFunctionBegin;
980   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
981   if (!ptap->P_oth) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
982 
983   ierr = MatZeroEntries(C);CHKERRQ(ierr);
984 
985   /* Get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
986   /*-----------------------------------------------------*/
987   if (ptap->reuse == MAT_REUSE_MATRIX) {
988     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
989     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
990   }
991 
992   po = (Mat_SeqAIJ*) p->B->data;
993   pd = (Mat_SeqAIJ*) p->A->data;
994 
995   ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr);
996   ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr);
997 
998   ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr);
999   ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr);
1000   ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr);
1001   cmaxr = 0;
1002   for (i=0; i<pon; i++) {
1003     cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]);
1004   }
1005   cd = (Mat_SeqAIJ*)(c->A)->data;
1006   co = (Mat_SeqAIJ*)(c->B)->data;
1007   for (i=0; i<pn; i++) {
1008     cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i]));
1009   }
1010   ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr);
1011   ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr);
1012   ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr);
1013   for (i=0; i<am && (pon || pn); i++) {
1014     ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr);
1015     nzi  = po->i[i+1] - po->i[i];
1016     dnzi = pd->i[i+1] - pd->i[i];
1017     if (!nzi && !dnzi) continue;
1018     ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr);
1019     voff = 0;
1020     ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr);
1021     if (!voff) continue;
1022 
1023     /* Form remote C(ii, :) */
1024     poj = po->j + po->i[i];
1025     poa = po->a + po->i[i];
1026     for (j=0; j<nzi; j++) {
1027       c_rmtjj = c_rmtj + ptap->c_rmti[poj[j]];
1028       c_rmtaa = c_rmta + ptap->c_rmti[poj[j]];
1029       for (jj=0; jj<voff; jj++) {
1030         apvaluestmp[jj] = apvalues[jj]*poa[j];
1031         /*If the row is empty */
1032         if (!c_rmtc[poj[j]]) {
1033           c_rmtjj[jj] = apindices[jj];
1034           c_rmtaa[jj] = apvaluestmp[jj];
1035           c_rmtc[poj[j]]++;
1036         } else {
1037           ierr = PetscFindInt(apindices[jj],c_rmtc[poj[j]],c_rmtjj,&loc);CHKERRQ(ierr);
1038           if (loc>=0){ /* hit */
1039             c_rmtaa[loc] += apvaluestmp[jj];
1040             ierr = PetscLogFlops(1.0);CHKERRQ(ierr);
1041           } else { /* new element */
1042             loc = -(loc+1);
1043             /* Move data backward */
1044             for (kk=c_rmtc[poj[j]]; kk>loc; kk--) {
1045               c_rmtjj[kk] = c_rmtjj[kk-1];
1046               c_rmtaa[kk] = c_rmtaa[kk-1];
1047             }/* End kk */
1048             c_rmtjj[loc] = apindices[jj];
1049             c_rmtaa[loc] = apvaluestmp[jj];
1050             c_rmtc[poj[j]]++;
1051           }
1052         }
1053       } /* End jj */
1054       ierr = PetscLogFlops(voff);CHKERRQ(ierr);
1055     } /* End j */
1056 
1057     /* Form local C(ii, :) */
1058     pdj = pd->j + pd->i[i];
1059     pda = pd->a + pd->i[i];
1060     for (j=0; j<dnzi; j++) {
1061       row = pcstart + pdj[j];
1062       for (jj=0; jj<voff; jj++) {
1063         apvaluestmp[jj] = apvalues[jj]*pda[j];
1064       }/* End kk */
1065       ierr = PetscLogFlops(voff);CHKERRQ(ierr);
1066       ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr);
1067     }/* End j */
1068   } /* End i */
1069 
1070   ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr);
1071   ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr);
1072   ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr);
1073 
1074   ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1075   ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr);
1076   ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_SUM);CHKERRQ(ierr);
1077   ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr);
1078   ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr);
1079 
1080   /* Add contributions from remote */
1081   for (i = 0; i < pn; i++) {
1082     row = i + pcstart;
1083     ierr = MatSetValues(C,1,&row,ptap->c_othi[i+1]-ptap->c_othi[i],c_othj+ptap->c_othi[i],c_otha+ptap->c_othi[i],ADD_VALUES);CHKERRQ(ierr);
1084   }
1085   ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr);
1086 
1087   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1088   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1089 
1090   ptap->reuse = MAT_REUSE_MATRIX;
1091 
1092   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
1093   if (ptap->freestruct) {
1094     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
1095   }
1096   PetscFunctionReturn(0);
1097 }
1098 
1099 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,PetscReal fill,Mat *C)
1100 {
1101   Mat_APMPI           *ptap;
1102   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*c;
1103   MPI_Comm            comm;
1104   Mat                 Cmpi;
1105   Mat_SeqAIJ          *pd,*po;
1106   MatType             mtype;
1107   PetscSF             sf;
1108   PetscSFNode         *iremote;
1109   PetscInt            rootspacesize,*rootspace,*rootspaceoffsets,nleaves;
1110   const PetscInt      *rootdegrees;
1111   PetscHSetI          ht,oht,*hta,*hto;
1112   PetscInt            pn,pon,*c_rmtc,i,j,nzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets;
1113   PetscInt            owner,lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj;
1114   PetscInt            nalg=2,alg=0;
1115   PetscBool           flg;
1116   const char          *algTypes[2] = {"overlapping","merged"};
1117   PetscErrorCode      ierr;
1118 
1119   PetscFunctionBegin;
1120   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1121 
1122   /* Create symbolic parallel matrix Cmpi */
1123   ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr);
1124   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1125   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1126   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
1127   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1128 
1129   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1130   ptap->reuse = MAT_INITIAL_MATRIX;
1131   ptap->algType = 2;
1132 
1133   /* Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1134   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1135 
1136   po = (Mat_SeqAIJ*)p->B->data;
1137   pd = (Mat_SeqAIJ*)p->A->data;
1138 
1139   /* This equals to the number of offdiag columns in P */
1140   ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr);
1141   /* offsets */
1142   ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr);
1143   /* The number of columns we will send to remote ranks */
1144   ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr);
1145   ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr);
1146   for (i=0; i<pon; i++) {
1147     ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr);
1148   }
1149   ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr);
1150   ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr);
1151   /* Create hash table to merge all columns for C(i, :) */
1152   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
1153 
1154   ptap->c_rmti[0] = 0;
1155   /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors)  */
1156   for (i=0; i<am && pon; i++) {
1157     /* Form one row of AP */
1158     ierr = PetscHSetIClear(ht);CHKERRQ(ierr);
1159     /* If the off diag is empty, we should not do any calculation */
1160     nzi = po->i[i+1] - po->i[i];
1161     if (!nzi) continue;
1162 
1163     ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,ht);CHKERRQ(ierr);
1164     ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr);
1165     /* If AP is empty, just continue */
1166     if (!htsize) continue;
1167     /* Form C(ii, :) */
1168     poj = po->j + po->i[i];
1169     for (j=0; j<nzi; j++) {
1170       ierr = PetscHSetIUpdate(hta[poj[j]],ht);CHKERRQ(ierr);
1171     }
1172   }
1173 
1174   for (i=0; i<pon; i++) {
1175     ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr);
1176     ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize;
1177     c_rmtc[i] = htsize;
1178   }
1179 
1180   ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr);
1181 
1182   for (i=0; i<pon; i++) {
1183     off = 0;
1184     ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr);
1185     ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr);
1186   }
1187   ierr = PetscFree(hta);CHKERRQ(ierr);
1188 
1189   ierr = PetscCalloc1(pon,&iremote);CHKERRQ(ierr);
1190   for (i=0; i<pon; i++) {
1191     owner = 0; lidx = 0;
1192     ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr);
1193     iremote[i].index = lidx;
1194     iremote[i].rank  = owner;
1195   }
1196 
1197   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
1198   ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1199   /* Reorder ranks properly so that the data handled by gather and scatter have the same order */
1200   ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr);
1201   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
1202   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1203   /* How many neighbors have contributions to my rows? */
1204   ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr);
1205   ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr);
1206   rootspacesize = 0;
1207   for (i = 0; i < pn; i++) {
1208     rootspacesize += rootdegrees[i];
1209   }
1210   ierr = PetscCalloc1(rootspacesize,&rootspace);CHKERRQ(ierr);
1211   ierr = PetscCalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr);
1212   /* Get information from leaves
1213    * Number of columns other people contribute to my rows
1214    * */
1215   ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr);
1216   ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr);
1217   ierr = PetscFree(c_rmtc);CHKERRQ(ierr);
1218   ierr = PetscCalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr);
1219   /* The number of columns is received for each row */
1220   ptap->c_othi[0] = 0;
1221   rootspacesize = 0;
1222   rootspaceoffsets[0] = 0;
1223   for (i = 0; i < pn; i++) {
1224     rcvncols = 0;
1225     for (j = 0; j<rootdegrees[i]; j++) {
1226       rcvncols += rootspace[rootspacesize];
1227       rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize];
1228       rootspacesize++;
1229     }
1230     ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols;
1231   }
1232   ierr = PetscFree(rootspace);CHKERRQ(ierr);
1233 
1234   ierr = PetscCalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr);
1235   ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr);
1236   ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr);
1237   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1238   ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr);
1239 
1240   ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr);
1241   nleaves = 0;
1242   for (i = 0; i<pon; i++) {
1243     owner = 0; lidx = 0;
1244     ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr);
1245     sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i];
1246     for (j=0; j<sendncols; j++) {
1247       iremote[nleaves].rank = owner;
1248       iremote[nleaves++].index = c_rmtoffsets[i] + j;
1249     }
1250   }
1251   ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr);
1252   ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr);
1253 
1254   ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr);
1255   ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1256   ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr);
1257   /* One to one map */
1258   ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1259 
1260   ierr = PetscCalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr);
1261   ierr = PetscHSetICreate(&oht);CHKERRQ(ierr);
1262   ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr);
1263   ierr = PetscMalloc2(pn,&hta,pn,&hto);CHKERRQ(ierr);
1264   for (i=0; i<pn; i++) {
1265     ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr);
1266     ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr);
1267   }
1268   /* Work on local part */
1269   /* 4) Pass 1: Estimate memory for C_loc */
1270   for (i=0; i<am && pn; i++) {
1271     ierr = PetscHSetIClear(ht);CHKERRQ(ierr);
1272     ierr = PetscHSetIClear(oht);CHKERRQ(ierr);
1273     nzi = pd->i[i+1] - pd->i[i];
1274     if (!nzi) continue;
1275 
1276     ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,oht);CHKERRQ(ierr);
1277     ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr);
1278     ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr);
1279     if (!(htsize+htosize)) continue;
1280     /* Form C(ii, :) */
1281     pdj = pd->j + pd->i[i];
1282     for (j=0; j<nzi; j++) {
1283       ierr = PetscHSetIUpdate(hta[pdj[j]],ht);CHKERRQ(ierr);
1284       ierr = PetscHSetIUpdate(hto[pdj[j]],oht);CHKERRQ(ierr);
1285     }
1286   }
1287 
1288   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
1289   ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr);
1290 
1291   /* Get remote data */
1292   ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1293   ierr = PetscFree(c_rmtj);CHKERRQ(ierr);
1294 
1295   for (i = 0; i < pn; i++) {
1296     nzi = ptap->c_othi[i+1] - ptap->c_othi[i];
1297     rdj = c_othj + ptap->c_othi[i];
1298     for (j = 0; j < nzi; j++) {
1299       col =  rdj[j];
1300       /* diag part */
1301       if (col>=pcstart && col<pcend) {
1302         ierr = PetscHSetIAdd(hta[i],col);CHKERRQ(ierr);
1303       } else { /* off diag */
1304         ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr);
1305       }
1306     }
1307     ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr);
1308     dnz[i] = htsize;
1309     ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr);
1310     ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr);
1311     onz[i] = htsize;
1312     ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr);
1313   }
1314 
1315   ierr = PetscFree2(hta,hto);CHKERRQ(ierr);
1316   ierr = PetscFree(c_othj);CHKERRQ(ierr);
1317 
1318   /* local sizes and preallocation */
1319   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1320   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
1321   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1322   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
1323 
1324   /* attach the supporting struct to Cmpi for reuse */
1325   c = (Mat_MPIAIJ*)Cmpi->data;
1326   c->ap           = ptap;
1327   ptap->duplicate = Cmpi->ops->duplicate;
1328   ptap->destroy   = Cmpi->ops->destroy;
1329   ptap->view      = Cmpi->ops->view;
1330 
1331   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1332   Cmpi->assembled        = PETSC_FALSE;
1333   /* pick an algorithm */
1334   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
1335   alg = 0;
1336   ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
1337   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1338   switch (alg) {
1339     case 0:
1340       Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce;
1341       break;
1342     case 1:
1343       Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged;
1344       break;
1345     default:
1346       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n");
1347   }
1348   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1349   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
1350   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1351   *C                     = Cmpi;
1352   PetscFunctionReturn(0);
1353 }
1354 
1355 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,PetscReal fill,Mat *C)
1356 {
1357   Mat_APMPI           *ptap;
1358   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*c;
1359   MPI_Comm            comm;
1360   Mat                 Cmpi;
1361   Mat_SeqAIJ          *pd,*po;
1362   MatType             mtype;
1363   PetscSF             sf;
1364   PetscSFNode         *iremote;
1365   PetscInt            rootspacesize,*rootspace,*rootspaceoffsets,nleaves;
1366   const PetscInt      *rootdegrees;
1367   PetscHSetI          ht,oht,*hta,*hto,*htd;
1368   PetscInt            pn,pon,*c_rmtc,i,j,nzi,dnzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets;
1369   PetscInt            owner,lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj;
1370   PetscInt            nalg=2,alg=0;
1371   PetscBool           flg;
1372   const char          *algTypes[2] = {"merged","overlapping"};
1373   PetscErrorCode      ierr;
1374 
1375   PetscFunctionBegin;
1376   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1377 
1378   /* Create symbolic parallel matrix Cmpi */
1379   ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr);
1380   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1381   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1382   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
1383   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1384 
1385   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1386   ptap->reuse = MAT_INITIAL_MATRIX;
1387   ptap->algType = 3;
1388 
1389   /* 0) Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1390   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1391 
1392   po = (Mat_SeqAIJ*)p->B->data;
1393   pd = (Mat_SeqAIJ*)p->A->data;
1394 
1395   /* This equals to the number of offdiag columns in P */
1396   ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr);
1397   /* offsets */
1398   ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr);
1399   /* The number of columns we will send to remote ranks */
1400   ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr);
1401   ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr);
1402   for (i=0; i<pon; i++) {
1403     ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr);
1404   }
1405   ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr);
1406   ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr);
1407   /* Create hash table to merge all columns for C(i, :) */
1408   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
1409   ierr = PetscHSetICreate(&oht);CHKERRQ(ierr);
1410   ierr = PetscMalloc2(pn,&htd,pn,&hto);CHKERRQ(ierr);
1411   for (i=0; i<pn; i++) {
1412     ierr = PetscHSetICreate(&htd[i]);CHKERRQ(ierr);
1413     ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr);
1414   }
1415   ptap->c_rmti[0] = 0;
1416   /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors)  */
1417   for (i=0; i<am && (pon || pn); i++) {
1418     /* Form one row of AP */
1419     ierr = PetscHSetIClear(ht);CHKERRQ(ierr);
1420     ierr = PetscHSetIClear(oht);CHKERRQ(ierr);
1421     /* If the off diag is empty, we should not do any calculation */
1422     nzi = po->i[i+1] - po->i[i];
1423     dnzi = pd->i[i+1] - pd->i[i];
1424     if (!nzi && !dnzi) continue;
1425 
1426     ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,oht);CHKERRQ(ierr);
1427     ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr);
1428     ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr);
1429     /* If AP is empty, just continue */
1430     if (!(htsize+htosize)) continue;
1431 
1432     /* Form remote C(ii, :) */
1433     poj = po->j + po->i[i];
1434     for (j=0; j<nzi; j++) {
1435       ierr = PetscHSetIUpdate(hta[poj[j]],ht);CHKERRQ(ierr);
1436       ierr = PetscHSetIUpdate(hta[poj[j]],oht);CHKERRQ(ierr);
1437     }
1438 
1439     /* Form local C(ii, :) */
1440     pdj = pd->j + pd->i[i];
1441     for (j=0; j<dnzi; j++) {
1442       ierr = PetscHSetIUpdate(htd[pdj[j]],ht);CHKERRQ(ierr);
1443       ierr = PetscHSetIUpdate(hto[pdj[j]],oht);CHKERRQ(ierr);
1444     }
1445   }
1446 
1447   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
1448   ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr);
1449 
1450   for (i=0; i<pon; i++) {
1451     ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr);
1452     ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize;
1453     c_rmtc[i] = htsize;
1454   }
1455 
1456   ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr);
1457 
1458   for (i=0; i<pon; i++) {
1459     off = 0;
1460     ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr);
1461     ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr);
1462   }
1463   ierr = PetscFree(hta);CHKERRQ(ierr);
1464 
1465   ierr = PetscCalloc1(pon,&iremote);CHKERRQ(ierr);
1466   for (i=0; i<pon; i++) {
1467     owner = 0; lidx = 0;
1468     ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr);
1469     iremote[i].index = lidx;
1470     iremote[i].rank  = owner;
1471   }
1472 
1473   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
1474   ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1475   /* Reorder ranks properly so that the data handled by gather and scatter have the same order */
1476   ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr);
1477   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
1478   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1479   /* How many neighbors have contributions to my rows? */
1480   ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr);
1481   ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr);
1482   rootspacesize = 0;
1483   for (i = 0; i < pn; i++) {
1484     rootspacesize += rootdegrees[i];
1485   }
1486   ierr = PetscCalloc1(rootspacesize,&rootspace);CHKERRQ(ierr);
1487   ierr = PetscCalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr);
1488   /* Get information from leaves
1489    * Number of columns other people contribute to my rows
1490    * */
1491   ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr);
1492   ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr);
1493   ierr = PetscFree(c_rmtc);CHKERRQ(ierr);
1494   ierr = PetscCalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr);
1495   /* The number of columns is received for each row */
1496   ptap->c_othi[0] = 0;
1497   rootspacesize = 0;
1498   rootspaceoffsets[0] = 0;
1499   for (i = 0; i < pn; i++) {
1500     rcvncols = 0;
1501     for (j = 0; j<rootdegrees[i]; j++) {
1502       rcvncols += rootspace[rootspacesize];
1503       rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize];
1504       rootspacesize++;
1505     }
1506     ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols;
1507   }
1508   ierr = PetscFree(rootspace);CHKERRQ(ierr);
1509 
1510   ierr = PetscCalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr);
1511   ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr);
1512   ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr);
1513   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1514   ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr);
1515 
1516   ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr);
1517   nleaves = 0;
1518   for (i = 0; i<pon; i++) {
1519     owner = 0; lidx = 0;
1520     ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr);
1521     sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i];
1522     for (j=0; j<sendncols; j++) {
1523       iremote[nleaves].rank = owner;
1524       iremote[nleaves++].index = c_rmtoffsets[i] + j;
1525     }
1526   }
1527   ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr);
1528   ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr);
1529 
1530   ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr);
1531   ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1532   ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr);
1533   /* One to one map */
1534   ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1535   /* Get remote data */
1536   ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1537   ierr = PetscFree(c_rmtj);CHKERRQ(ierr);
1538   ierr = PetscCalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr);
1539   ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr);
1540 
1541   for (i = 0; i < pn; i++) {
1542     nzi = ptap->c_othi[i+1] - ptap->c_othi[i];
1543     rdj = c_othj + ptap->c_othi[i];
1544     for (j = 0; j < nzi; j++) {
1545       col =  rdj[j];
1546       /* diag part */
1547       if (col>=pcstart && col<pcend) {
1548         ierr = PetscHSetIAdd(htd[i],col);CHKERRQ(ierr);
1549       } else { /* off diag */
1550         ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr);
1551       }
1552     }
1553     ierr = PetscHSetIGetSize(htd[i],&htsize);CHKERRQ(ierr);
1554     dnz[i] = htsize;
1555     ierr = PetscHSetIDestroy(&htd[i]);CHKERRQ(ierr);
1556     ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr);
1557     onz[i] = htsize;
1558     ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr);
1559   }
1560 
1561   ierr = PetscFree2(htd,hto);CHKERRQ(ierr);
1562   ierr = PetscFree(c_othj);CHKERRQ(ierr);
1563 
1564   /* local sizes and preallocation */
1565   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1566   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
1567   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1568   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
1569 
1570   /* attach the supporting struct to Cmpi for reuse */
1571   c = (Mat_MPIAIJ*)Cmpi->data;
1572   c->ap           = ptap;
1573   ptap->duplicate = Cmpi->ops->duplicate;
1574   ptap->destroy   = Cmpi->ops->destroy;
1575   ptap->view      = Cmpi->ops->view;
1576 
1577   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1578   Cmpi->assembled        = PETSC_FALSE;
1579   /* pick an algorithm */
1580   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
1581   alg = 0;
1582   ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
1583   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1584   switch (alg) {
1585     case 0:
1586       Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged;
1587       break;
1588     case 1:
1589       Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce;
1590       break;
1591     default:
1592       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n");
1593   }
1594   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1595   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
1596   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1597   *C                     = Cmpi;
1598   PetscFunctionReturn(0);
1599 }
1600 
1601 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
1602 {
1603   PetscErrorCode      ierr;
1604   Mat_APMPI           *ptap;
1605   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
1606   MPI_Comm            comm;
1607   PetscMPIInt         size,rank;
1608   Mat                 Cmpi;
1609   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1610   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
1611   PetscInt            *lnk,i,k,pnz,row,nsend;
1612   PetscBT             lnkbt;
1613   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
1614   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1615   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
1616   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1617   MPI_Request         *swaits,*rwaits;
1618   MPI_Status          *sstatus,rstatus;
1619   PetscLayout         rowmap;
1620   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
1621   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
1622   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
1623   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
1624   PetscScalar         *apv;
1625   PetscTable          ta;
1626   MatType             mtype;
1627   const char          *prefix;
1628 #if defined(PETSC_USE_INFO)
1629   PetscReal           apfill;
1630 #endif
1631 
1632   PetscFunctionBegin;
1633   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1634   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1635   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1636 
1637   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
1638 
1639   /* create symbolic parallel matrix Cmpi */
1640   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1641   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1642   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
1643 
1644   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
1645   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
1646 
1647   /* create struct Mat_APMPI and attached it to C later */
1648   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1649   ptap->reuse = MAT_INITIAL_MATRIX;
1650   ptap->algType = 1;
1651 
1652   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1653   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1654   /* get P_loc by taking all local rows of P */
1655   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1656 
1657   /* (0) compute Rd = Pd^T, Ro = Po^T  */
1658   /* --------------------------------- */
1659   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1660   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1661 
1662   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
1663   /* ----------------------------------------------------------------- */
1664   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1665   if (ptap->P_oth) p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1666 
1667   /* create and initialize a linked list */
1668   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
1669   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
1670   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
1671   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
1672   /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */
1673 
1674   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
1675 
1676   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
1677   if (ao) {
1678     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
1679   } else {
1680     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
1681   }
1682   current_space = free_space;
1683   nspacedouble  = 0;
1684 
1685   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
1686   api[0] = 0;
1687   for (i=0; i<am; i++) {
1688     /* diagonal portion: Ad[i,:]*P */
1689     ai = ad->i; pi = p_loc->i;
1690     nzi = ai[i+1] - ai[i];
1691     aj  = ad->j + ai[i];
1692     for (j=0; j<nzi; j++) {
1693       row  = aj[j];
1694       pnz  = pi[row+1] - pi[row];
1695       Jptr = p_loc->j + pi[row];
1696       /* add non-zero cols of P into the sorted linked list lnk */
1697       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1698     }
1699     /* off-diagonal portion: Ao[i,:]*P */
1700     if (ao) {
1701       ai = ao->i; pi = p_oth->i;
1702       nzi = ai[i+1] - ai[i];
1703       aj  = ao->j + ai[i];
1704       for (j=0; j<nzi; j++) {
1705         row  = aj[j];
1706         pnz  = pi[row+1] - pi[row];
1707         Jptr = p_oth->j + pi[row];
1708         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1709       }
1710     }
1711     apnz     = lnk[0];
1712     api[i+1] = api[i] + apnz;
1713     if (ap_rmax < apnz) ap_rmax = apnz;
1714 
1715     /* if free space is not available, double the total space in the list */
1716     if (current_space->local_remaining<apnz) {
1717       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1718       nspacedouble++;
1719     }
1720 
1721     /* Copy data into free space, then initialize lnk */
1722     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1723 
1724     current_space->array           += apnz;
1725     current_space->local_used      += apnz;
1726     current_space->local_remaining -= apnz;
1727   }
1728   /* Allocate space for apj and apv, initialize apj, and */
1729   /* destroy list of free space and other temporary array(s) */
1730   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
1731   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
1732   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1733 
1734   /* Create AP_loc for reuse */
1735   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
1736 
1737 #if defined(PETSC_USE_INFO)
1738   if (ao) {
1739     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
1740   } else {
1741     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
1742   }
1743   ptap->AP_loc->info.mallocs           = nspacedouble;
1744   ptap->AP_loc->info.fill_ratio_given  = fill;
1745   ptap->AP_loc->info.fill_ratio_needed = apfill;
1746 
1747   if (api[am]) {
1748     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);
1749     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
1750   } else {
1751     ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
1752   }
1753 #endif
1754 
1755   /* (2-1) compute symbolic Co = Ro*AP_loc  */
1756   /* ------------------------------------ */
1757   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1758   ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr);
1759   ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr);
1760   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
1761 
1762   /* (3) send coj of C_oth to other processors  */
1763   /* ------------------------------------------ */
1764   /* determine row ownership */
1765   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
1766   rowmap->n  = pn;
1767   rowmap->bs = 1;
1768   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
1769   owners = rowmap->range;
1770 
1771   /* determine the number of messages to send, their lengths */
1772   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
1773   /* TODO: Use PetscCalloc2() */
1774   ierr = PetscArrayzero(len_s,size);CHKERRQ(ierr);
1775   ierr = PetscArrayzero(len_si,size);CHKERRQ(ierr);
1776 
1777   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1778   coi   = c_oth->i; coj = c_oth->j;
1779   con   = ptap->C_oth->rmap->n;
1780   proc  = 0;
1781   for (i=0; i<con; i++) {
1782     while (prmap[i] >= owners[proc+1]) proc++;
1783     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
1784     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1785   }
1786 
1787   len          = 0; /* max length of buf_si[], see (4) */
1788   owners_co[0] = 0;
1789   nsend        = 0;
1790   for (proc=0; proc<size; proc++) {
1791     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1792     if (len_s[proc]) {
1793       nsend++;
1794       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1795       len         += len_si[proc];
1796     }
1797   }
1798 
1799   /* determine the number and length of messages to receive for coi and coj  */
1800   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
1801   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
1802 
1803   /* post the Irecv and Isend of coj */
1804   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1805   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1806   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
1807   for (proc=0, k=0; proc<size; proc++) {
1808     if (!len_s[proc]) continue;
1809     i    = owners_co[proc];
1810     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1811     k++;
1812   }
1813 
1814   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
1815   /* ---------------------------------------- */
1816   ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr);
1817   ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr);
1818   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
1819   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
1820 
1821   /* receives coj are complete */
1822   for (i=0; i<nrecv; i++) {
1823     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1824   }
1825   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1826   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
1827 
1828   /* add received column indices into ta to update Crmax */
1829   for (k=0; k<nrecv; k++) {/* k-th received message */
1830     Jptr = buf_rj[k];
1831     for (j=0; j<len_r[k]; j++) {
1832       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
1833     }
1834   }
1835   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
1836   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
1837 
1838   /* (4) send and recv coi */
1839   /*-----------------------*/
1840   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1841   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1842   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1843   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1844   for (proc=0,k=0; proc<size; proc++) {
1845     if (!len_s[proc]) continue;
1846     /* form outgoing message for i-structure:
1847          buf_si[0]:                 nrows to be sent
1848                [1:nrows]:           row index (global)
1849                [nrows+1:2*nrows+1]: i-structure index
1850     */
1851     /*-------------------------------------------*/
1852     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1853     buf_si_i    = buf_si + nrows+1;
1854     buf_si[0]   = nrows;
1855     buf_si_i[0] = 0;
1856     nrows       = 0;
1857     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1858       nzi = coi[i+1] - coi[i];
1859       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1860       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1861       nrows++;
1862     }
1863     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1864     k++;
1865     buf_si += len_si[proc];
1866   }
1867   for (i=0; i<nrecv; i++) {
1868     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1869   }
1870   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1871   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
1872 
1873   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
1874   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1875   ierr = PetscFree(swaits);CHKERRQ(ierr);
1876   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1877 
1878   /* (5) compute the local portion of Cmpi      */
1879   /* ------------------------------------------ */
1880   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
1881   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
1882   current_space = free_space;
1883 
1884   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
1885   for (k=0; k<nrecv; k++) {
1886     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1887     nrows       = *buf_ri_k[k];
1888     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1889     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1890   }
1891 
1892   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
1893   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
1894   for (i=0; i<pn; i++) {
1895     /* add C_loc into Cmpi */
1896     nzi  = c_loc->i[i+1] - c_loc->i[i];
1897     Jptr = c_loc->j + c_loc->i[i];
1898     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1899 
1900     /* add received col data into lnk */
1901     for (k=0; k<nrecv; k++) { /* k-th received message */
1902       if (i == *nextrow[k]) { /* i-th row */
1903         nzi  = *(nextci[k]+1) - *nextci[k];
1904         Jptr = buf_rj[k] + *nextci[k];
1905         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1906         nextrow[k]++; nextci[k]++;
1907       }
1908     }
1909     nzi = lnk[0];
1910 
1911     /* copy data into free space, then initialize lnk */
1912     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1913     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
1914   }
1915   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1916   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1917   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
1918 
1919   /* local sizes and preallocation */
1920   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1921   if (P->cmap->bs > 0) {
1922     ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr);
1923     ierr = PetscLayoutSetBlockSize(Cmpi->cmap,P->cmap->bs);CHKERRQ(ierr);
1924   }
1925   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1926   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1927 
1928   /* members in merge */
1929   ierr = PetscFree(id_r);CHKERRQ(ierr);
1930   ierr = PetscFree(len_r);CHKERRQ(ierr);
1931   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
1932   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
1933   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
1934   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
1935   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
1936 
1937   /* attach the supporting struct to Cmpi for reuse */
1938   c = (Mat_MPIAIJ*)Cmpi->data;
1939   c->ap           = ptap;
1940   ptap->duplicate = Cmpi->ops->duplicate;
1941   ptap->destroy   = Cmpi->ops->destroy;
1942   ptap->view      = Cmpi->ops->view;
1943   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
1944 
1945   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1946   Cmpi->assembled        = PETSC_FALSE;
1947   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1948   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
1949   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1950   *C                     = Cmpi;
1951   PetscFunctionReturn(0);
1952 }
1953 
1954 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
1955 {
1956   PetscErrorCode    ierr;
1957   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1958   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1959   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
1960   Mat_APMPI         *ptap = c->ap;
1961   Mat               AP_loc,C_loc,C_oth;
1962   PetscInt          i,rstart,rend,cm,ncols,row;
1963   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
1964   PetscScalar       *apa;
1965   const PetscInt    *cols;
1966   const PetscScalar *vals;
1967 
1968   PetscFunctionBegin;
1969   if (!ptap->AP_loc) {
1970     MPI_Comm comm;
1971     ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1972     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1973   }
1974 
1975   ierr = MatZeroEntries(C);CHKERRQ(ierr);
1976   /* 1) get R = Pd^T,Ro = Po^T */
1977   if (ptap->reuse == MAT_REUSE_MATRIX) {
1978     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1979     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1980   }
1981 
1982   /* 2) get AP_loc */
1983   AP_loc = ptap->AP_loc;
1984   ap = (Mat_SeqAIJ*)AP_loc->data;
1985 
1986   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
1987   /*-----------------------------------------------------*/
1988   if (ptap->reuse == MAT_REUSE_MATRIX) {
1989     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
1990     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1991     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1992   }
1993 
1994   /* 2-2) compute numeric A_loc*P - dominating part */
1995   /* ---------------------------------------------- */
1996   /* get data from symbolic products */
1997   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1998   if (ptap->P_oth) {
1999     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
2000   }
2001   apa   = ptap->apa;
2002   api   = ap->i;
2003   apj   = ap->j;
2004   for (i=0; i<am; i++) {
2005     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
2006     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
2007     apnz = api[i+1] - api[i];
2008     for (j=0; j<apnz; j++) {
2009       col = apj[j+api[i]];
2010       ap->a[j+ap->i[i]] = apa[col];
2011       apa[col] = 0.0;
2012     }
2013   }
2014 
2015   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
2016   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
2017   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
2018   C_loc = ptap->C_loc;
2019   C_oth = ptap->C_oth;
2020 
2021   /* add C_loc and Co to to C */
2022   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
2023 
2024   /* C_loc -> C */
2025   cm    = C_loc->rmap->N;
2026   c_seq = (Mat_SeqAIJ*)C_loc->data;
2027   cols = c_seq->j;
2028   vals = c_seq->a;
2029 
2030 
2031   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
2032   /* when there are no off-processor parts.  */
2033   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
2034   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
2035   /* a table, and other, more complex stuff has to be done. */
2036   if (C->assembled) {
2037     C->was_assembled = PETSC_TRUE;
2038     C->assembled     = PETSC_FALSE;
2039   }
2040   if (C->was_assembled) {
2041     for (i=0; i<cm; i++) {
2042       ncols = c_seq->i[i+1] - c_seq->i[i];
2043       row = rstart + i;
2044       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
2045       cols += ncols; vals += ncols;
2046     }
2047   } else {
2048     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
2049   }
2050 
2051   /* Co -> C, off-processor part */
2052   cm = C_oth->rmap->N;
2053   c_seq = (Mat_SeqAIJ*)C_oth->data;
2054   cols = c_seq->j;
2055   vals = c_seq->a;
2056   for (i=0; i<cm; i++) {
2057     ncols = c_seq->i[i+1] - c_seq->i[i];
2058     row = p->garray[i];
2059     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
2060     cols += ncols; vals += ncols;
2061   }
2062 
2063   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2064   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2065 
2066   ptap->reuse = MAT_REUSE_MATRIX;
2067 
2068   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
2069   if (ptap->freestruct) {
2070     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
2071   }
2072   PetscFunctionReturn(0);
2073 }
2074