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