xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 63c55180f3a1699b45d3541b604931f10ad6dae3)
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 = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr);
270     AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa);
271     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
272   }
273   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,api[AP_loc->rmap->n],apj,&nout,apj);CHKERRQ(ierr);
274   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);
275 
276   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
277   /* Always use scalable version since we are in the MPI scalable version */
278   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
279   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
280 
281   C_loc = ptap->C_loc;
282   C_oth = ptap->C_oth;
283 
284   /* add C_loc and Co to to C */
285   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
286 
287   /* C_loc -> C */
288   cm    = C_loc->rmap->N;
289   c_seq = (Mat_SeqAIJ*)C_loc->data;
290   cols = c_seq->j;
291   vals = c_seq->a;
292   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_loc->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr);
293 
294   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
295   /* when there are no off-processor parts.  */
296   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
297   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
298   /* a table, and other, more complex stuff has to be done. */
299   if (C->assembled) {
300     C->was_assembled = PETSC_TRUE;
301     C->assembled     = PETSC_FALSE;
302   }
303   if (C->was_assembled) {
304     for (i=0; i<cm; i++) {
305       ncols = c_seq->i[i+1] - c_seq->i[i];
306       row = rstart + i;
307       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
308       cols += ncols; vals += ncols;
309     }
310   } else {
311     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
312   }
313   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_loc->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr);
314   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);
315 
316   /* Co -> C, off-processor part */
317   cm = C_oth->rmap->N;
318   c_seq = (Mat_SeqAIJ*)C_oth->data;
319   cols = c_seq->j;
320   vals = c_seq->a;
321   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_oth->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr);
322   for (i=0; i<cm; i++) {
323     ncols = c_seq->i[i+1] - c_seq->i[i];
324     row = p->garray[i];
325     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
326     cols += ncols; vals += ncols;
327   }
328   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
329   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
330 
331   ptap->reuse = MAT_REUSE_MATRIX;
332 
333   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_oth->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr);
334   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);
335 
336   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
337   if (ptap->freestruct) {
338     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
339   }
340   PetscFunctionReturn(0);
341 }
342 
343 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
344 {
345   PetscErrorCode      ierr;
346   Mat_APMPI           *ptap;
347   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
348   MPI_Comm            comm;
349   PetscMPIInt         size,rank;
350   Mat                 Cmpi,P_loc,P_oth;
351   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
352   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
353   PetscInt            *lnk,i,k,pnz,row,nsend;
354   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
355   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
356   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
357   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
358   MPI_Request         *swaits,*rwaits;
359   MPI_Status          *sstatus,rstatus;
360   PetscLayout         rowmap;
361   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
362   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
363   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi,nout;
364   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
365   PetscScalar         *apv;
366   PetscTable          ta;
367   MatType             mtype;
368   const char          *prefix;
369 #if defined(PETSC_USE_INFO)
370   PetscReal           apfill;
371 #endif
372 
373   PetscFunctionBegin;
374   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
375   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
376   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
377 
378   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
379 
380   /* create symbolic parallel matrix Cmpi */
381   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
382   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
383   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
384 
385   /* create struct Mat_APMPI and attached it to C later */
386   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
387   ptap->reuse = MAT_INITIAL_MATRIX;
388   ptap->algType = 0;
389 
390   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
391   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
392   /* get P_loc by taking all local rows of P */
393   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
394 
395   ptap->P_loc = P_loc;
396   ptap->P_oth = P_oth;
397 
398   /* (0) compute Rd = Pd^T, Ro = Po^T  */
399   /* --------------------------------- */
400   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
401   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
402 
403   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
404   /* ----------------------------------------------------------------- */
405   p_loc  = (Mat_SeqAIJ*)P_loc->data;
406   if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data;
407 
408   /* create and initialize a linked list */
409   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
410   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
411   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
412   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
413 
414   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
415 
416   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
417   if (ao) {
418     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
419   } else {
420     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
421   }
422   current_space = free_space;
423   nspacedouble  = 0;
424 
425   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
426   api[0] = 0;
427   for (i=0; i<am; i++) {
428     /* diagonal portion: Ad[i,:]*P */
429     ai = ad->i; pi = p_loc->i;
430     nzi = ai[i+1] - ai[i];
431     aj  = ad->j + ai[i];
432     for (j=0; j<nzi; j++) {
433       row  = aj[j];
434       pnz  = pi[row+1] - pi[row];
435       Jptr = p_loc->j + pi[row];
436       /* add non-zero cols of P into the sorted linked list lnk */
437       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
438     }
439     /* off-diagonal portion: Ao[i,:]*P */
440     if (ao) {
441       ai = ao->i; pi = p_oth->i;
442       nzi = ai[i+1] - ai[i];
443       aj  = ao->j + ai[i];
444       for (j=0; j<nzi; j++) {
445         row  = aj[j];
446         pnz  = pi[row+1] - pi[row];
447         Jptr = p_oth->j + pi[row];
448         ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
449       }
450     }
451     apnz     = lnk[0];
452     api[i+1] = api[i] + apnz;
453 
454     /* if free space is not available, double the total space in the list */
455     if (current_space->local_remaining<apnz) {
456       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
457       nspacedouble++;
458     }
459 
460     /* Copy data into free space, then initialize lnk */
461     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
462 
463     current_space->array           += apnz;
464     current_space->local_used      += apnz;
465     current_space->local_remaining -= apnz;
466   }
467   /* Allocate space for apj and apv, initialize apj, and */
468   /* destroy list of free space and other temporary array(s) */
469   ierr = PetscCalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
470   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
471   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
472 
473   /* Create AP_loc for reuse */
474   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
475   ierr = MatSeqAIJCompactOutExtraColumns_SeqAIJ(ptap->AP_loc, &ptap->ltog);CHKERRQ(ierr);
476 
477 #if defined(PETSC_USE_INFO)
478   if (ao) {
479     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
480   } else {
481     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
482   }
483   ptap->AP_loc->info.mallocs           = nspacedouble;
484   ptap->AP_loc->info.fill_ratio_given  = fill;
485   ptap->AP_loc->info.fill_ratio_needed = apfill;
486 
487   if (api[am]) {
488     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);
489     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
490   } else {
491     ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
492   }
493 #endif
494 
495   /* (2-1) compute symbolic Co = Ro*AP_loc  */
496   /* ------------------------------------ */
497   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
498   ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr);
499   ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr);
500   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
501 
502   /* (3) send coj of C_oth to other processors  */
503   /* ------------------------------------------ */
504   /* determine row ownership */
505   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
506   rowmap->n  = pn;
507   rowmap->bs = 1;
508   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
509   owners = rowmap->range;
510 
511   /* determine the number of messages to send, their lengths */
512   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
513   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
514   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));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   }
784   for (j=0; j<nzi; j++) {
785     ra   = aa[j];
786     row  = aj[j];
787     nzpi = po->i[row+1] - po->i[row];
788     pj = po->j + po->i[row];
789     pa = po->a + po->i[row];
790     for (k=0; k<nzpi; k++) {
791       ierr = PetscHMapIVAddValue(hmap,p->garray[pj[k]],ra*pa[k]);CHKERRQ(ierr);
792     }
793   }
794 
795 
796   /* off diagonal part: Ao[i, :]*P_oth */
797   if (ao) {
798     ai = ao->i;
799     pi = p_oth->i;
800     nzi = ai[i+1] - ai[i];
801     aj  = ao->j + ai[i];
802     aa  = ao->a + ai[i];
803     for (j=0; j<nzi; j++) {
804       row  = aj[j];
805       ra   = aa[j];
806       pnz  = pi[row+1] - pi[row];
807       p_othcols = p_oth->j + pi[row];
808       pa   = p_oth->a + pi[row];
809       for (col=0; col<pnz; col++) {
810         ierr = PetscHMapIVAddValue(hmap,p_othcols[col],ra*pa[col]);CHKERRQ(ierr);
811       }
812     }
813   } /* end if (ao) */
814   PetscFunctionReturn(0);
815 }
816 
817 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,Mat C)
818 {
819   PetscErrorCode    ierr;
820   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
821   Mat_SeqAIJ        *cd,*co,*po,*pd;
822   Mat_APMPI         *ptap = c->ap;
823   PetscHMapIV       hmap;
824   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;
825   PetscScalar       *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa;
826   MPI_Comm          comm;
827 
828   PetscFunctionBegin;
829   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
830   if (!ptap->P_oth) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
831 
832   ierr = MatZeroEntries(C);CHKERRQ(ierr);
833 
834   /* Get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
835   /*-----------------------------------------------------*/
836   if (ptap->reuse == MAT_REUSE_MATRIX) {
837     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
838     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
839   }
840 
841   po = (Mat_SeqAIJ*) p->B->data;
842   pd = (Mat_SeqAIJ*) p->A->data;
843 
844   ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr);
845 
846   ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr);
847   ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr);
848   cmaxr = 0;
849   for (i=0; i<pon; i++) {
850     cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]);
851   }
852   ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr);
853   ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr);
854   ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr);
855   for (i=0; i<am && pon; i++) {
856     ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr);
857     nzi = po->i[i+1] - po->i[i];
858     if (!nzi) continue;
859     ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr);
860     voff = 0;
861     ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr);
862     if (!voff) continue;
863     /*ierr = PetscMemzero(c_rmtc,sizeof(PetscInt)*pon);CHKERRQ(ierr);*/
864 
865     /* Form C(ii, :) */
866     poj = po->j + po->i[i];
867     poa = po->a + po->i[i];
868     for (j=0; j<nzi; j++) {
869       c_rmtjj = c_rmtj + ptap->c_rmti[poj[j]];
870       c_rmtaa = c_rmta + ptap->c_rmti[poj[j]];
871       for (jj=0; jj<voff; jj++) {
872         apvaluestmp[jj] = apvalues[jj]*poa[j];
873         /*If the row is empty */
874         if (!c_rmtc[poj[j]]) {
875           c_rmtjj[jj] = apindices[jj];
876           c_rmtaa[jj] = apvaluestmp[jj];
877           c_rmtc[poj[j]]++;
878         } else {
879           ierr = PetscFindInt(apindices[jj],c_rmtc[poj[j]],c_rmtjj,&loc);CHKERRQ(ierr);
880           if (loc>=0){ /* hit */
881             c_rmtaa[loc] += apvaluestmp[jj];
882           } else { /* new element */
883             loc = -(loc+1);
884             /* Move data backward */
885             for (kk=c_rmtc[poj[j]]; kk>loc; kk--) {
886               c_rmtjj[kk] = c_rmtjj[kk-1];
887               c_rmtaa[kk] = c_rmtaa[kk-1];
888             }/* End kk */
889             c_rmtjj[loc] = apindices[jj];
890             c_rmtaa[loc] = apvaluestmp[jj];
891             c_rmtc[poj[j]]++;
892           }
893         }
894       } /* End jj */
895     } /* End j */
896   } /* End i */
897 
898   ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr);
899 
900   ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr);
901   ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr);
902 
903   ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
904   ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr);
905   ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr);
906   cd = (Mat_SeqAIJ*)(c->A)->data;
907   co = (Mat_SeqAIJ*)(c->B)->data;
908 
909   cmaxr = 0;
910   for (i=0; i<pn; i++) {
911     cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i]));
912   }
913   ierr = PetscCalloc5(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pn,&dcc,pn,&occ);CHKERRQ(ierr);
914   ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr);
915   ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr);
916   for (i=0; i<am && pn; i++) {
917     ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr);
918     nzi = pd->i[i+1] - pd->i[i];
919     if (!nzi) continue;
920     ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr);
921     voff = 0;
922     ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr);
923     if (!voff) continue;
924     /* Form C(ii, :) */
925     pdj = pd->j + pd->i[i];
926     pda = pd->a + pd->i[i];
927     for (j=0; j<nzi; j++) {
928       row = pcstart + pdj[j];
929       for (jj=0; jj<voff; jj++) {
930         apvaluestmp[jj] = apvalues[jj]*pda[j];
931       }
932       ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr);
933     }
934   }
935 
936   ierr = MatGetOwnershipRangeColumn(C,&ccstart,&ccend);CHKERRQ(ierr);
937   ierr = PetscFree5(apindices,apvalues,apvaluestmp,dcc,occ);CHKERRQ(ierr);
938   ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr);
939   ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_SUM);CHKERRQ(ierr);
940   ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr);
941   ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr);
942 
943   /* Add contributions from remote */
944   for (i = 0; i < pn; i++) {
945     row = i + pcstart;
946     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);
947   }
948   ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr);
949 
950   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
951   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
952 
953   ptap->reuse = MAT_REUSE_MATRIX;
954 
955   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
956   if (ptap->freestruct) {
957     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
958   }
959   PetscFunctionReturn(0);
960 }
961 
962 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,Mat C)
963 {
964   PetscErrorCode    ierr;
965   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
966   Mat_SeqAIJ        *cd,*co,*po,*pd;
967   Mat_APMPI         *ptap = c->ap;
968   PetscHMapIV       hmap;
969   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;
970   PetscScalar       *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa;
971   MPI_Comm          comm;
972 
973   PetscFunctionBegin;
974   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
975   if (!ptap->P_oth) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
976 
977   ierr = MatZeroEntries(C);CHKERRQ(ierr);
978 
979   /* Get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
980   /*-----------------------------------------------------*/
981   if (ptap->reuse == MAT_REUSE_MATRIX) {
982     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
983     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
984   }
985 
986   po = (Mat_SeqAIJ*) p->B->data;
987   pd = (Mat_SeqAIJ*) p->A->data;
988 
989   ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr);
990   ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr);
991 
992   ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr);
993   ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr);
994   ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr);
995   cmaxr = 0;
996   for (i=0; i<pon; i++) {
997     cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]);
998   }
999   cd = (Mat_SeqAIJ*)(c->A)->data;
1000   co = (Mat_SeqAIJ*)(c->B)->data;
1001   for (i=0; i<pn; i++) {
1002     cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i]));
1003   }
1004   ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr);
1005   ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr);
1006   ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr);
1007   for (i=0; i<am && (pon || pn); i++) {
1008     ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr);
1009     nzi  = po->i[i+1] - po->i[i];
1010     dnzi = pd->i[i+1] - pd->i[i];
1011     if (!nzi && !dnzi) continue;
1012     ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr);
1013     voff = 0;
1014     ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr);
1015     if (!voff) continue;
1016 
1017     /* Form remote C(ii, :) */
1018     poj = po->j + po->i[i];
1019     poa = po->a + po->i[i];
1020     for (j=0; j<nzi; j++) {
1021       c_rmtjj = c_rmtj + ptap->c_rmti[poj[j]];
1022       c_rmtaa = c_rmta + ptap->c_rmti[poj[j]];
1023       for (jj=0; jj<voff; jj++) {
1024         apvaluestmp[jj] = apvalues[jj]*poa[j];
1025         /*If the row is empty */
1026         if (!c_rmtc[poj[j]]) {
1027           c_rmtjj[jj] = apindices[jj];
1028           c_rmtaa[jj] = apvaluestmp[jj];
1029           c_rmtc[poj[j]]++;
1030         } else {
1031           ierr = PetscFindInt(apindices[jj],c_rmtc[poj[j]],c_rmtjj,&loc);CHKERRQ(ierr);
1032           if (loc>=0){ /* hit */
1033             c_rmtaa[loc] += apvaluestmp[jj];
1034           } else { /* new element */
1035             loc = -(loc+1);
1036             /* Move data backward */
1037             for (kk=c_rmtc[poj[j]]; kk>loc; kk--) {
1038               c_rmtjj[kk] = c_rmtjj[kk-1];
1039               c_rmtaa[kk] = c_rmtaa[kk-1];
1040             }/* End kk */
1041             c_rmtjj[loc] = apindices[jj];
1042             c_rmtaa[loc] = apvaluestmp[jj];
1043             c_rmtc[poj[j]]++;
1044           }
1045         }
1046       } /* End jj */
1047     } /* End j */
1048 
1049     /* Form local C(ii, :) */
1050     pdj = pd->j + pd->i[i];
1051     pda = pd->a + pd->i[i];
1052     for (j=0; j<dnzi; j++) {
1053       row = pcstart + pdj[j];
1054       for (jj=0; jj<voff; jj++) {
1055         apvaluestmp[jj] = apvalues[jj]*pda[j];
1056       }/* End kk */
1057       ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr);
1058     }/* End j */
1059   } /* End i */
1060 
1061   ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr);
1062   ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr);
1063   ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr);
1064 
1065   ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1066   ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr);
1067   ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_SUM);CHKERRQ(ierr);
1068   ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr);
1069   ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr);
1070 
1071   /* Add contributions from remote */
1072   for (i = 0; i < pn; i++) {
1073     row = i + pcstart;
1074     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);
1075   }
1076   ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr);
1077 
1078   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1079   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1080 
1081   ptap->reuse = MAT_REUSE_MATRIX;
1082 
1083   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
1084   if (ptap->freestruct) {
1085     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
1086   }
1087   PetscFunctionReturn(0);
1088 }
1089 
1090 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,PetscReal fill,Mat *C)
1091 {
1092   Mat_APMPI           *ptap;
1093   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*c;
1094   MPI_Comm            comm;
1095   Mat                 Cmpi;
1096   Mat_SeqAIJ          *pd,*po;
1097   MatType             mtype;
1098   PetscSF             sf;
1099   PetscSFNode         *iremote;
1100   PetscInt            rootspacesize,*rootspace,*rootspaceoffsets,nleaves;
1101   const PetscInt      *rootdegrees;
1102   PetscHSetI          ht,oht,*hta,*hto;
1103   PetscInt            pn,pon,*c_rmtc,i,j,nzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets;
1104   PetscInt            owner,lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj;
1105   PetscInt            nalg=2,alg=0;
1106   PetscBool           flg;
1107   const char          *algTypes[2] = {"overlapping","merged"};
1108   PetscErrorCode      ierr;
1109 
1110   PetscFunctionBegin;
1111   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1112 
1113   /* Create symbolic parallel matrix Cmpi */
1114   ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr);
1115   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1116   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1117   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
1118   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1119 
1120   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1121   ptap->reuse = MAT_INITIAL_MATRIX;
1122   ptap->algType = 2;
1123 
1124   /* Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1125   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1126 
1127   po = (Mat_SeqAIJ*)p->B->data;
1128   pd = (Mat_SeqAIJ*)p->A->data;
1129 
1130   /* This equals to the number of offdiag columns in P */
1131   ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr);
1132   /* offsets */
1133   ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr);
1134   /* The number of columns we will send to remote ranks */
1135   ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr);
1136   ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr);
1137   for (i=0; i<pon; i++) {
1138     ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr);
1139   }
1140   ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr);
1141   ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr);
1142   /* Create hash table to merge all columns for C(i, :) */
1143   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
1144 
1145   ptap->c_rmti[0] = 0;
1146   /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors)  */
1147   for (i=0; i<am && pon; i++) {
1148     /* Form one row of AP */
1149     ierr = PetscHSetIClear(ht);CHKERRQ(ierr);
1150     /* If the off diag is empty, we should not do any calculation */
1151     nzi = po->i[i+1] - po->i[i];
1152     if (!nzi) continue;
1153 
1154     ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,ht);CHKERRQ(ierr);
1155     ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr);
1156     /* If AP is empty, just continue */
1157     if (!htsize) continue;
1158     /* Form C(ii, :) */
1159     poj = po->j + po->i[i];
1160     for (j=0; j<nzi; j++) {
1161       ierr = PetscHSetIUpdate(hta[poj[j]],ht);CHKERRQ(ierr);
1162     }
1163   }
1164 
1165   for (i=0; i<pon; i++) {
1166     ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr);
1167     ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize;
1168     c_rmtc[i] = htsize;
1169   }
1170 
1171   ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr);
1172 
1173   for (i=0; i<pon; i++) {
1174     off = 0;
1175     ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr);
1176     ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr);
1177   }
1178   ierr = PetscFree(hta);CHKERRQ(ierr);
1179 
1180   ierr = PetscCalloc1(pon,&iremote);CHKERRQ(ierr);
1181   for (i=0; i<pon; i++) {
1182     owner = 0; lidx = 0;
1183     ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr);
1184     iremote[i].index = lidx;
1185     iremote[i].rank  = owner;
1186   }
1187 
1188   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
1189   ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1190   /* Reorder ranks properly so that the data handled by gather and scatter have the same order */
1191   ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr);
1192   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
1193   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1194   /* How many neighbors have contributions to my rows? */
1195   ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr);
1196   ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr);
1197   rootspacesize = 0;
1198   for (i = 0; i < pn; i++) {
1199     rootspacesize += rootdegrees[i];
1200   }
1201   ierr = PetscCalloc1(rootspacesize,&rootspace);CHKERRQ(ierr);
1202   ierr = PetscCalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr);
1203   /* Get information from leaves
1204    * Number of columns other people contribute to my rows
1205    * */
1206   ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr);
1207   ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr);
1208   ierr = PetscFree(c_rmtc);CHKERRQ(ierr);
1209   ierr = PetscCalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr);
1210   /* The number of columns is received for each row */
1211   ptap->c_othi[0] = 0;
1212   rootspacesize = 0;
1213   rootspaceoffsets[0] = 0;
1214   for (i = 0; i < pn; i++) {
1215     rcvncols = 0;
1216     for (j = 0; j<rootdegrees[i]; j++) {
1217       rcvncols += rootspace[rootspacesize];
1218       rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize];
1219       rootspacesize++;
1220     }
1221     ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols;
1222   }
1223   ierr = PetscFree(rootspace);CHKERRQ(ierr);
1224 
1225   ierr = PetscCalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr);
1226   ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr);
1227   ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr);
1228   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1229   ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr);
1230 
1231   ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr);
1232   nleaves = 0;
1233   for (i = 0; i<pon; i++) {
1234     owner = 0; lidx = 0;
1235     ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr);
1236     sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i];
1237     for (j=0; j<sendncols; j++) {
1238       iremote[nleaves].rank = owner;
1239       iremote[nleaves++].index = c_rmtoffsets[i] + j;
1240     }
1241   }
1242   ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr);
1243   ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr);
1244 
1245   ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr);
1246   ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1247   ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr);
1248   /* One to one map */
1249   ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1250 
1251   ierr = PetscCalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr);
1252   ierr = PetscHSetICreate(&oht);CHKERRQ(ierr);
1253   ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr);
1254   ierr = PetscMalloc2(pn,&hta,pn,&hto);CHKERRQ(ierr);
1255   for (i=0; i<pn; i++) {
1256     ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr);
1257     ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr);
1258   }
1259   /* Work on local part */
1260   /* 4) Pass 1: Estimate memory for C_loc */
1261   for (i=0; i<am && pn; i++) {
1262     ierr = PetscHSetIClear(ht);CHKERRQ(ierr);
1263     ierr = PetscHSetIClear(oht);CHKERRQ(ierr);
1264     nzi = pd->i[i+1] - pd->i[i];
1265     if (!nzi) continue;
1266 
1267     ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,oht);CHKERRQ(ierr);
1268     ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr);
1269     ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr);
1270     if (!(htsize+htosize)) continue;
1271     /* Form C(ii, :) */
1272     pdj = pd->j + pd->i[i];
1273     for (j=0; j<nzi; j++) {
1274       ierr = PetscHSetIUpdate(hta[pdj[j]],ht);CHKERRQ(ierr);
1275       ierr = PetscHSetIUpdate(hto[pdj[j]],oht);CHKERRQ(ierr);
1276     }
1277   }
1278 
1279   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
1280   ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr);
1281 
1282   /* Get remote data */
1283   ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1284   ierr = PetscFree(c_rmtj);CHKERRQ(ierr);
1285 
1286   for (i = 0; i < pn; i++) {
1287     nzi = ptap->c_othi[i+1] - ptap->c_othi[i];
1288     rdj = c_othj + ptap->c_othi[i];
1289     for (j = 0; j < nzi; j++) {
1290       col =  rdj[j];
1291       /* diag part */
1292       if (col>=pcstart && col<pcend) {
1293         ierr = PetscHSetIAdd(hta[i],col);CHKERRQ(ierr);
1294       } else { /* off diag */
1295         ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr);
1296       }
1297     }
1298     ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr);
1299     dnz[i] = htsize;
1300     ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr);
1301     ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr);
1302     onz[i] = htsize;
1303     ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr);
1304   }
1305 
1306   ierr = PetscFree2(hta,hto);CHKERRQ(ierr);
1307   ierr = PetscFree(c_othj);CHKERRQ(ierr);
1308 
1309   /* local sizes and preallocation */
1310   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1311   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
1312   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1313   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
1314 
1315   /* attach the supporting struct to Cmpi for reuse */
1316   c = (Mat_MPIAIJ*)Cmpi->data;
1317   c->ap           = ptap;
1318   ptap->duplicate = Cmpi->ops->duplicate;
1319   ptap->destroy   = Cmpi->ops->destroy;
1320   ptap->view      = Cmpi->ops->view;
1321 
1322   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1323   Cmpi->assembled        = PETSC_FALSE;
1324   /* pick an algorithm */
1325   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
1326   alg = 0;
1327   ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
1328   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1329   switch (alg) {
1330     case 0:
1331       Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce;
1332       break;
1333     case 1:
1334       Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged;
1335       break;
1336     default:
1337       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n");
1338   }
1339   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1340   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
1341   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1342   *C                     = Cmpi;
1343   PetscFunctionReturn(0);
1344 }
1345 
1346 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,PetscReal fill,Mat *C)
1347 {
1348   Mat_APMPI           *ptap;
1349   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*c;
1350   MPI_Comm            comm;
1351   Mat                 Cmpi;
1352   Mat_SeqAIJ          *pd,*po;
1353   MatType             mtype;
1354   PetscSF             sf;
1355   PetscSFNode         *iremote;
1356   PetscInt            rootspacesize,*rootspace,*rootspaceoffsets,nleaves;
1357   const PetscInt      *rootdegrees;
1358   PetscHSetI          ht,oht,*hta,*hto,*htd;
1359   PetscInt            pn,pon,*c_rmtc,i,j,nzi,dnzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets;
1360   PetscInt            owner,lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj;
1361   PetscInt            nalg=2,alg=0;
1362   PetscBool           flg;
1363   const char          *algTypes[2] = {"merged","overlapping"};
1364   PetscErrorCode      ierr;
1365 
1366   PetscFunctionBegin;
1367   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1368 
1369   /* Create symbolic parallel matrix Cmpi */
1370   ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr);
1371   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1372   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1373   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
1374   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1375 
1376   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1377   ptap->reuse = MAT_INITIAL_MATRIX;
1378   ptap->algType = 3;
1379 
1380   /* 0) Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1381   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1382 
1383   po = (Mat_SeqAIJ*)p->B->data;
1384   pd = (Mat_SeqAIJ*)p->A->data;
1385 
1386   /* This equals to the number of offdiag columns in P */
1387   ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr);
1388   /* offsets */
1389   ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr);
1390   /* The number of columns we will send to remote ranks */
1391   ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr);
1392   ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr);
1393   for (i=0; i<pon; i++) {
1394     ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr);
1395   }
1396   ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr);
1397   ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr);
1398   /* Create hash table to merge all columns for C(i, :) */
1399   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
1400   ierr = PetscHSetICreate(&oht);CHKERRQ(ierr);
1401   ierr = PetscMalloc2(pn,&htd,pn,&hto);CHKERRQ(ierr);
1402   for (i=0; i<pn; i++) {
1403     ierr = PetscHSetICreate(&htd[i]);CHKERRQ(ierr);
1404     ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr);
1405   }
1406   ptap->c_rmti[0] = 0;
1407   /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors)  */
1408   for (i=0; i<am && (pon || pn); i++) {
1409     /* Form one row of AP */
1410     ierr = PetscHSetIClear(ht);CHKERRQ(ierr);
1411     ierr = PetscHSetIClear(oht);CHKERRQ(ierr);
1412     /* If the off diag is empty, we should not do any calculation */
1413     nzi = po->i[i+1] - po->i[i];
1414     dnzi = pd->i[i+1] - pd->i[i];
1415     if (!nzi && !dnzi) continue;
1416 
1417     ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,oht);CHKERRQ(ierr);
1418     ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr);
1419     ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr);
1420     /* If AP is empty, just continue */
1421     if (!(htsize+htosize)) continue;
1422 
1423     /* Form remote C(ii, :) */
1424     poj = po->j + po->i[i];
1425     for (j=0; j<nzi; j++) {
1426       ierr = PetscHSetIUpdate(hta[poj[j]],ht);CHKERRQ(ierr);
1427       ierr = PetscHSetIUpdate(hta[poj[j]],oht);CHKERRQ(ierr);
1428     }
1429 
1430     /* Form local C(ii, :) */
1431     pdj = pd->j + pd->i[i];
1432     for (j=0; j<dnzi; j++) {
1433       ierr = PetscHSetIUpdate(htd[pdj[j]],ht);CHKERRQ(ierr);
1434       ierr = PetscHSetIUpdate(hto[pdj[j]],oht);CHKERRQ(ierr);
1435     }
1436   }
1437 
1438   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
1439   ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr);
1440 
1441   for (i=0; i<pon; i++) {
1442     ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr);
1443     ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize;
1444     c_rmtc[i] = htsize;
1445   }
1446 
1447   ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr);
1448 
1449   for (i=0; i<pon; i++) {
1450     off = 0;
1451     ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr);
1452     ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr);
1453   }
1454   ierr = PetscFree(hta);CHKERRQ(ierr);
1455 
1456   ierr = PetscCalloc1(pon,&iremote);CHKERRQ(ierr);
1457   for (i=0; i<pon; i++) {
1458     owner = 0; lidx = 0;
1459     ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr);
1460     iremote[i].index = lidx;
1461     iremote[i].rank  = owner;
1462   }
1463 
1464   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
1465   ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1466   /* Reorder ranks properly so that the data handled by gather and scatter have the same order */
1467   ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr);
1468   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
1469   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1470   /* How many neighbors have contributions to my rows? */
1471   ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr);
1472   ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr);
1473   rootspacesize = 0;
1474   for (i = 0; i < pn; i++) {
1475     rootspacesize += rootdegrees[i];
1476   }
1477   ierr = PetscCalloc1(rootspacesize,&rootspace);CHKERRQ(ierr);
1478   ierr = PetscCalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr);
1479   /* Get information from leaves
1480    * Number of columns other people contribute to my rows
1481    * */
1482   ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr);
1483   ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr);
1484   ierr = PetscFree(c_rmtc);CHKERRQ(ierr);
1485   ierr = PetscCalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr);
1486   /* The number of columns is received for each row */
1487   ptap->c_othi[0] = 0;
1488   rootspacesize = 0;
1489   rootspaceoffsets[0] = 0;
1490   for (i = 0; i < pn; i++) {
1491     rcvncols = 0;
1492     for (j = 0; j<rootdegrees[i]; j++) {
1493       rcvncols += rootspace[rootspacesize];
1494       rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize];
1495       rootspacesize++;
1496     }
1497     ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols;
1498   }
1499   ierr = PetscFree(rootspace);CHKERRQ(ierr);
1500 
1501   ierr = PetscCalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr);
1502   ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr);
1503   ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr);
1504   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1505   ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr);
1506 
1507   ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr);
1508   nleaves = 0;
1509   for (i = 0; i<pon; i++) {
1510     owner = 0; lidx = 0;
1511     ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr);
1512     sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i];
1513     for (j=0; j<sendncols; j++) {
1514       iremote[nleaves].rank = owner;
1515       iremote[nleaves++].index = c_rmtoffsets[i] + j;
1516     }
1517   }
1518   ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr);
1519   ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr);
1520 
1521   ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr);
1522   ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1523   ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr);
1524   /* One to one map */
1525   ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1526   /* Get remote data */
1527   ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr);
1528   ierr = PetscFree(c_rmtj);CHKERRQ(ierr);
1529   ierr = PetscCalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr);
1530   ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr);
1531 
1532   for (i = 0; i < pn; i++) {
1533     nzi = ptap->c_othi[i+1] - ptap->c_othi[i];
1534     rdj = c_othj + ptap->c_othi[i];
1535     for (j = 0; j < nzi; j++) {
1536       col =  rdj[j];
1537       /* diag part */
1538       if (col>=pcstart && col<pcend) {
1539         ierr = PetscHSetIAdd(htd[i],col);CHKERRQ(ierr);
1540       } else { /* off diag */
1541         ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr);
1542       }
1543     }
1544     ierr = PetscHSetIGetSize(htd[i],&htsize);CHKERRQ(ierr);
1545     dnz[i] = htsize;
1546     ierr = PetscHSetIDestroy(&htd[i]);CHKERRQ(ierr);
1547     ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr);
1548     onz[i] = htsize;
1549     ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr);
1550   }
1551 
1552   ierr = PetscFree2(htd,hto);CHKERRQ(ierr);
1553   ierr = PetscFree(c_othj);CHKERRQ(ierr);
1554 
1555   /* local sizes and preallocation */
1556   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1557   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
1558   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1559   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
1560 
1561   /* attach the supporting struct to Cmpi for reuse */
1562   c = (Mat_MPIAIJ*)Cmpi->data;
1563   c->ap           = ptap;
1564   ptap->duplicate = Cmpi->ops->duplicate;
1565   ptap->destroy   = Cmpi->ops->destroy;
1566   ptap->view      = Cmpi->ops->view;
1567 
1568   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1569   Cmpi->assembled        = PETSC_FALSE;
1570   /* pick an algorithm */
1571   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
1572   alg = 0;
1573   ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
1574   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1575   switch (alg) {
1576     case 0:
1577       Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged;
1578       break;
1579     case 1:
1580       Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce;
1581       break;
1582     default:
1583       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n");
1584   }
1585   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1586   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
1587   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1588   *C                     = Cmpi;
1589   PetscFunctionReturn(0);
1590 }
1591 
1592 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
1593 {
1594   PetscErrorCode      ierr;
1595   Mat_APMPI           *ptap;
1596   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
1597   MPI_Comm            comm;
1598   PetscMPIInt         size,rank;
1599   Mat                 Cmpi;
1600   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1601   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
1602   PetscInt            *lnk,i,k,pnz,row,nsend;
1603   PetscBT             lnkbt;
1604   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
1605   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1606   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
1607   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1608   MPI_Request         *swaits,*rwaits;
1609   MPI_Status          *sstatus,rstatus;
1610   PetscLayout         rowmap;
1611   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
1612   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
1613   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
1614   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
1615   PetscScalar         *apv;
1616   PetscTable          ta;
1617   MatType             mtype;
1618   const char          *prefix;
1619 #if defined(PETSC_USE_INFO)
1620   PetscReal           apfill;
1621 #endif
1622 
1623   PetscFunctionBegin;
1624   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1625   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1626   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1627 
1628   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
1629 
1630   /* create symbolic parallel matrix Cmpi */
1631   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1632   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1633   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
1634 
1635   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
1636   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
1637 
1638   /* create struct Mat_APMPI and attached it to C later */
1639   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1640   ptap->reuse = MAT_INITIAL_MATRIX;
1641   ptap->algType = 1;
1642 
1643   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1644   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1645   /* get P_loc by taking all local rows of P */
1646   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1647 
1648   /* (0) compute Rd = Pd^T, Ro = Po^T  */
1649   /* --------------------------------- */
1650   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1651   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1652 
1653   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
1654   /* ----------------------------------------------------------------- */
1655   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1656   if (ptap->P_oth) p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1657 
1658   /* create and initialize a linked list */
1659   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
1660   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
1661   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
1662   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
1663   /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */
1664 
1665   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
1666 
1667   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
1668   if (ao) {
1669     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
1670   } else {
1671     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
1672   }
1673   current_space = free_space;
1674   nspacedouble  = 0;
1675 
1676   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
1677   api[0] = 0;
1678   for (i=0; i<am; i++) {
1679     /* diagonal portion: Ad[i,:]*P */
1680     ai = ad->i; pi = p_loc->i;
1681     nzi = ai[i+1] - ai[i];
1682     aj  = ad->j + ai[i];
1683     for (j=0; j<nzi; j++) {
1684       row  = aj[j];
1685       pnz  = pi[row+1] - pi[row];
1686       Jptr = p_loc->j + pi[row];
1687       /* add non-zero cols of P into the sorted linked list lnk */
1688       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1689     }
1690     /* off-diagonal portion: Ao[i,:]*P */
1691     if (ao) {
1692       ai = ao->i; pi = p_oth->i;
1693       nzi = ai[i+1] - ai[i];
1694       aj  = ao->j + ai[i];
1695       for (j=0; j<nzi; j++) {
1696         row  = aj[j];
1697         pnz  = pi[row+1] - pi[row];
1698         Jptr = p_oth->j + pi[row];
1699         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1700       }
1701     }
1702     apnz     = lnk[0];
1703     api[i+1] = api[i] + apnz;
1704     if (ap_rmax < apnz) ap_rmax = apnz;
1705 
1706     /* if free space is not available, double the total space in the list */
1707     if (current_space->local_remaining<apnz) {
1708       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1709       nspacedouble++;
1710     }
1711 
1712     /* Copy data into free space, then initialize lnk */
1713     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1714 
1715     current_space->array           += apnz;
1716     current_space->local_used      += apnz;
1717     current_space->local_remaining -= apnz;
1718   }
1719   /* Allocate space for apj and apv, initialize apj, and */
1720   /* destroy list of free space and other temporary array(s) */
1721   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
1722   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
1723   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1724 
1725   /* Create AP_loc for reuse */
1726   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
1727 
1728 #if defined(PETSC_USE_INFO)
1729   if (ao) {
1730     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
1731   } else {
1732     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
1733   }
1734   ptap->AP_loc->info.mallocs           = nspacedouble;
1735   ptap->AP_loc->info.fill_ratio_given  = fill;
1736   ptap->AP_loc->info.fill_ratio_needed = apfill;
1737 
1738   if (api[am]) {
1739     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);
1740     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
1741   } else {
1742     ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
1743   }
1744 #endif
1745 
1746   /* (2-1) compute symbolic Co = Ro*AP_loc  */
1747   /* ------------------------------------ */
1748   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1749   ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr);
1750   ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr);
1751   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
1752 
1753   /* (3) send coj of C_oth to other processors  */
1754   /* ------------------------------------------ */
1755   /* determine row ownership */
1756   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
1757   rowmap->n  = pn;
1758   rowmap->bs = 1;
1759   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
1760   owners = rowmap->range;
1761 
1762   /* determine the number of messages to send, their lengths */
1763   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
1764   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1765   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1766 
1767   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1768   coi   = c_oth->i; coj = c_oth->j;
1769   con   = ptap->C_oth->rmap->n;
1770   proc  = 0;
1771   for (i=0; i<con; i++) {
1772     while (prmap[i] >= owners[proc+1]) proc++;
1773     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
1774     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1775   }
1776 
1777   len          = 0; /* max length of buf_si[], see (4) */
1778   owners_co[0] = 0;
1779   nsend        = 0;
1780   for (proc=0; proc<size; proc++) {
1781     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1782     if (len_s[proc]) {
1783       nsend++;
1784       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1785       len         += len_si[proc];
1786     }
1787   }
1788 
1789   /* determine the number and length of messages to receive for coi and coj  */
1790   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
1791   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
1792 
1793   /* post the Irecv and Isend of coj */
1794   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1795   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1796   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
1797   for (proc=0, k=0; proc<size; proc++) {
1798     if (!len_s[proc]) continue;
1799     i    = owners_co[proc];
1800     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1801     k++;
1802   }
1803 
1804   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
1805   /* ---------------------------------------- */
1806   ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr);
1807   ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr);
1808   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
1809   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
1810 
1811   /* receives coj are complete */
1812   for (i=0; i<nrecv; i++) {
1813     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1814   }
1815   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1816   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
1817 
1818   /* add received column indices into ta to update Crmax */
1819   for (k=0; k<nrecv; k++) {/* k-th received message */
1820     Jptr = buf_rj[k];
1821     for (j=0; j<len_r[k]; j++) {
1822       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
1823     }
1824   }
1825   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
1826   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
1827 
1828   /* (4) send and recv coi */
1829   /*-----------------------*/
1830   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1831   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1832   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1833   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1834   for (proc=0,k=0; proc<size; proc++) {
1835     if (!len_s[proc]) continue;
1836     /* form outgoing message for i-structure:
1837          buf_si[0]:                 nrows to be sent
1838                [1:nrows]:           row index (global)
1839                [nrows+1:2*nrows+1]: i-structure index
1840     */
1841     /*-------------------------------------------*/
1842     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1843     buf_si_i    = buf_si + nrows+1;
1844     buf_si[0]   = nrows;
1845     buf_si_i[0] = 0;
1846     nrows       = 0;
1847     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1848       nzi = coi[i+1] - coi[i];
1849       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1850       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1851       nrows++;
1852     }
1853     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1854     k++;
1855     buf_si += len_si[proc];
1856   }
1857   for (i=0; i<nrecv; i++) {
1858     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1859   }
1860   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1861   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
1862 
1863   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
1864   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1865   ierr = PetscFree(swaits);CHKERRQ(ierr);
1866   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1867 
1868   /* (5) compute the local portion of Cmpi      */
1869   /* ------------------------------------------ */
1870   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
1871   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
1872   current_space = free_space;
1873 
1874   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
1875   for (k=0; k<nrecv; k++) {
1876     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1877     nrows       = *buf_ri_k[k];
1878     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1879     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1880   }
1881 
1882   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
1883   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
1884   for (i=0; i<pn; i++) {
1885     /* add C_loc into Cmpi */
1886     nzi  = c_loc->i[i+1] - c_loc->i[i];
1887     Jptr = c_loc->j + c_loc->i[i];
1888     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1889 
1890     /* add received col data into lnk */
1891     for (k=0; k<nrecv; k++) { /* k-th received message */
1892       if (i == *nextrow[k]) { /* i-th row */
1893         nzi  = *(nextci[k]+1) - *nextci[k];
1894         Jptr = buf_rj[k] + *nextci[k];
1895         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1896         nextrow[k]++; nextci[k]++;
1897       }
1898     }
1899     nzi = lnk[0];
1900 
1901     /* copy data into free space, then initialize lnk */
1902     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1903     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
1904   }
1905   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1906   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1907   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
1908 
1909   /* local sizes and preallocation */
1910   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1911   if (P->cmap->bs > 0) {
1912     ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr);
1913     ierr = PetscLayoutSetBlockSize(Cmpi->cmap,P->cmap->bs);CHKERRQ(ierr);
1914   }
1915   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1916   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1917 
1918   /* members in merge */
1919   ierr = PetscFree(id_r);CHKERRQ(ierr);
1920   ierr = PetscFree(len_r);CHKERRQ(ierr);
1921   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
1922   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
1923   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
1924   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
1925   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
1926 
1927   /* attach the supporting struct to Cmpi for reuse */
1928   c = (Mat_MPIAIJ*)Cmpi->data;
1929   c->ap           = ptap;
1930   ptap->duplicate = Cmpi->ops->duplicate;
1931   ptap->destroy   = Cmpi->ops->destroy;
1932   ptap->view      = Cmpi->ops->view;
1933   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
1934 
1935   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1936   Cmpi->assembled        = PETSC_FALSE;
1937   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1938   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
1939   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1940   *C                     = Cmpi;
1941   PetscFunctionReturn(0);
1942 }
1943 
1944 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
1945 {
1946   PetscErrorCode    ierr;
1947   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1948   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1949   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
1950   Mat_APMPI         *ptap = c->ap;
1951   Mat               AP_loc,C_loc,C_oth;
1952   PetscInt          i,rstart,rend,cm,ncols,row;
1953   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
1954   PetscScalar       *apa;
1955   const PetscInt    *cols;
1956   const PetscScalar *vals;
1957 
1958   PetscFunctionBegin;
1959   if (!ptap->AP_loc) {
1960     MPI_Comm comm;
1961     ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1962     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1963   }
1964 
1965   ierr = MatZeroEntries(C);CHKERRQ(ierr);
1966   /* 1) get R = Pd^T,Ro = Po^T */
1967   if (ptap->reuse == MAT_REUSE_MATRIX) {
1968     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1969     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1970   }
1971 
1972   /* 2) get AP_loc */
1973   AP_loc = ptap->AP_loc;
1974   ap = (Mat_SeqAIJ*)AP_loc->data;
1975 
1976   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
1977   /*-----------------------------------------------------*/
1978   if (ptap->reuse == MAT_REUSE_MATRIX) {
1979     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
1980     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1981     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1982   }
1983 
1984   /* 2-2) compute numeric A_loc*P - dominating part */
1985   /* ---------------------------------------------- */
1986   /* get data from symbolic products */
1987   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1988   if (ptap->P_oth) {
1989     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1990   }
1991   apa   = ptap->apa;
1992   api   = ap->i;
1993   apj   = ap->j;
1994   for (i=0; i<am; i++) {
1995     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
1996     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
1997     apnz = api[i+1] - api[i];
1998     for (j=0; j<apnz; j++) {
1999       col = apj[j+api[i]];
2000       ap->a[j+ap->i[i]] = apa[col];
2001       apa[col] = 0.0;
2002     }
2003     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
2004   }
2005 
2006   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
2007   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
2008   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
2009   C_loc = ptap->C_loc;
2010   C_oth = ptap->C_oth;
2011 
2012   /* add C_loc and Co to to C */
2013   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
2014 
2015   /* C_loc -> C */
2016   cm    = C_loc->rmap->N;
2017   c_seq = (Mat_SeqAIJ*)C_loc->data;
2018   cols = c_seq->j;
2019   vals = c_seq->a;
2020 
2021 
2022   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
2023   /* when there are no off-processor parts.  */
2024   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
2025   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
2026   /* a table, and other, more complex stuff has to be done. */
2027   if (C->assembled) {
2028     C->was_assembled = PETSC_TRUE;
2029     C->assembled     = PETSC_FALSE;
2030   }
2031   if (C->was_assembled) {
2032     for (i=0; i<cm; i++) {
2033       ncols = c_seq->i[i+1] - c_seq->i[i];
2034       row = rstart + i;
2035       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
2036       cols += ncols; vals += ncols;
2037     }
2038   } else {
2039     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
2040   }
2041 
2042   /* Co -> C, off-processor part */
2043   cm = C_oth->rmap->N;
2044   c_seq = (Mat_SeqAIJ*)C_oth->data;
2045   cols = c_seq->j;
2046   vals = c_seq->a;
2047   for (i=0; i<cm; i++) {
2048     ncols = c_seq->i[i+1] - c_seq->i[i];
2049     row = p->garray[i];
2050     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
2051     cols += ncols; vals += ncols;
2052   }
2053 
2054   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2055   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2056 
2057   ptap->reuse = MAT_REUSE_MATRIX;
2058 
2059   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
2060   if (ptap->freestruct) {
2061     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
2062   }
2063   PetscFunctionReturn(0);
2064 }
2065