xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 7601faf09de2bc280184d9c73acff348d5eb2a25)
1 
2 /*
3   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
4           C = A * B
5 */
6 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
7 #include <../src/mat/utils/freespace.h>
8 #include <../src/mat/impls/aij/mpi/mpiaij.h>
9 #include <petscbt.h>
10 #include <../src/mat/impls/dense/mpi/mpidense.h>
11 #include <petsc-private/vecimpl.h>
12 
13 #undef __FUNCT__
14 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
15 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
16 {
17   PetscErrorCode ierr;
18 
19   PetscFunctionBegin;
20   if (scall == MAT_INITIAL_MATRIX) {
21     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
22   }
23 
24   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
25   PetscFunctionReturn(0);
26 }
27 
28 #undef __FUNCT__
29 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
30 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
31 {
32   PetscErrorCode ierr;
33   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
34   Mat_PtAPMPI    *ptap = a->ptap;
35 
36   PetscFunctionBegin;
37   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
38   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
39   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
40   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
41   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
42   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
43   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
44   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
45   ierr = ptap->destroy(A);CHKERRQ(ierr);
46   ierr = PetscFree(ptap);CHKERRQ(ierr);
47   PetscFunctionReturn(0);
48 }
49 
50 #undef __FUNCT__
51 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
52 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
53 {
54   PetscErrorCode ierr;
55   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
56   Mat_PtAPMPI    *ptap = a->ptap;
57 
58   PetscFunctionBegin;
59   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
60 
61   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
62   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
68 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
69 {
70   PetscErrorCode ierr;
71   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
72   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
73   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
74   PetscInt       *adi=ad->i,*adj,*aoi=ao->i,*aoj;
75   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
76   Mat_SeqAIJ     *p_loc,*p_oth;
77   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
78   PetscScalar    *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
79   PetscInt       cm   =C->rmap->n,anz,pnz;
80   Mat_PtAPMPI    *ptap=c->ptap;
81   PetscInt       *api,*apj,*apJ,i,j,k,row;
82   PetscInt       cstart=C->cmap->rstart;
83   PetscInt       cdnz,conz,k0,k1;
84 
85   PetscFunctionBegin;
86   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
87   /*-----------------------------------------------------*/
88   /* update numerical values of P_oth and P_loc */
89   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
90   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
91 
92   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
93   /*----------------------------------------------------------*/
94   /* get data from symbolic products */
95   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
96   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
97   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
98   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
99 
100   /* get apa for storing dense row A[i,:]*P */
101   apa = ptap->apa;
102 
103   api = ptap->api;
104   apj = ptap->apj;
105   for (i=0; i<cm; i++) {
106     /* diagonal portion of A */
107     anz = adi[i+1] - adi[i];
108     adj = ad->j + adi[i];
109     ada = ad->a + adi[i];
110     for (j=0; j<anz; j++) {
111       row = adj[j];
112       pnz = pi_loc[row+1] - pi_loc[row];
113       pj  = pj_loc + pi_loc[row];
114       pa  = pa_loc + pi_loc[row];
115 
116       /* perform dense axpy */
117       valtmp = ada[j];
118       for (k=0; k<pnz; k++) {
119         apa[pj[k]] += valtmp*pa[k];
120       }
121       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
122     }
123 
124     /* off-diagonal portion of A */
125     anz = aoi[i+1] - aoi[i];
126     aoj = ao->j + aoi[i];
127     aoa = ao->a + aoi[i];
128     for (j=0; j<anz; j++) {
129       row = aoj[j];
130       pnz = pi_oth[row+1] - pi_oth[row];
131       pj  = pj_oth + pi_oth[row];
132       pa  = pa_oth + pi_oth[row];
133 
134       /* perform dense axpy */
135       valtmp = aoa[j];
136       for (k=0; k<pnz; k++) {
137         apa[pj[k]] += valtmp*pa[k];
138       }
139       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
140     }
141 
142     /* set values in C */
143     apJ  = apj + api[i];
144     cdnz = cd->i[i+1] - cd->i[i];
145     conz = co->i[i+1] - co->i[i];
146 
147     /* 1st off-diagoanl part of C */
148     ca = coa + co->i[i];
149     k  = 0;
150     for (k0=0; k0<conz; k0++) {
151       if (apJ[k] >= cstart) break;
152       ca[k0]      = apa[apJ[k]];
153       apa[apJ[k]] = 0.0;
154       k++;
155     }
156 
157     /* diagonal part of C */
158     ca = cda + cd->i[i];
159     for (k1=0; k1<cdnz; k1++) {
160       ca[k1]      = apa[apJ[k]];
161       apa[apJ[k]] = 0.0;
162       k++;
163     }
164 
165     /* 2nd off-diagoanl part of C */
166     ca = coa + co->i[i];
167     for (; k0<conz; k0++) {
168       ca[k0]      = apa[apJ[k]];
169       apa[apJ[k]] = 0.0;
170       k++;
171     }
172   }
173   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
174   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
175   PetscFunctionReturn(0);
176 }
177 
178 #undef __FUNCT__
179 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
180 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
181 {
182   PetscErrorCode     ierr;
183   MPI_Comm           comm;
184   Mat                Cmpi;
185   Mat_PtAPMPI        *ptap;
186   PetscFreeSpaceList free_space=NULL,current_space=NULL;
187   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
188   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
189   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
190   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
191   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
192   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
193   PetscBT            lnkbt;
194   PetscScalar        *apa;
195   PetscReal          afill;
196   PetscBool          scalable=PETSC_TRUE;
197   PetscInt           nlnk_max,armax,prmax;
198 
199   PetscFunctionBegin;
200   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
201   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
202     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
203   }
204 
205   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
206   ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,NULL);CHKERRQ(ierr);
207   if (scalable) {
208     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);CHKERRQ(ierr);
209     PetscFunctionReturn(0);
210   }
211   ierr = PetscOptionsEnd();CHKERRQ(ierr);
212 
213   /* create struct Mat_PtAPMPI and attached it to C later */
214   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
215 
216   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
217   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
218 
219   /* get P_loc by taking all local rows of P */
220   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
221 
222   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
223   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
224   pi_loc = p_loc->i; pj_loc = p_loc->j;
225   pi_oth = p_oth->i; pj_oth = p_oth->j;
226 
227   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
228   /*-------------------------------------------------------------------*/
229   ierr      = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
230   ptap->api = api;
231   api[0]    = 0;
232 
233   /* create and initialize a linked list */
234   armax    = ad->rmax+ao->rmax;
235   prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
236   nlnk_max = armax*prmax;
237   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
238   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
239 
240   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
241   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
242 
243   current_space = free_space;
244 
245   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
246   for (i=0; i<am; i++) {
247     /* diagonal portion of A */
248     nzi = adi[i+1] - adi[i];
249     for (j=0; j<nzi; j++) {
250       row  = *adj++;
251       pnz  = pi_loc[row+1] - pi_loc[row];
252       Jptr = pj_loc + pi_loc[row];
253       /* add non-zero cols of P into the sorted linked list lnk */
254       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
255     }
256     /* off-diagonal portion of A */
257     nzi = aoi[i+1] - aoi[i];
258     for (j=0; j<nzi; j++) {
259       row  = *aoj++;
260       pnz  = pi_oth[row+1] - pi_oth[row];
261       Jptr = pj_oth + pi_oth[row];
262       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
263     }
264 
265     apnz     = lnk[0];
266     api[i+1] = api[i] + apnz;
267 
268     /* if free space is not available, double the total space in the list */
269     if (current_space->local_remaining<apnz) {
270       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
271       nspacedouble++;
272     }
273 
274     /* Copy data into free space, then initialize lnk */
275     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
276     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
277 
278     current_space->array           += apnz;
279     current_space->local_used      += apnz;
280     current_space->local_remaining -= apnz;
281   }
282 
283   /* Allocate space for apj, initialize apj, and */
284   /* destroy list of free space and other temporary array(s) */
285   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
286   apj  = ptap->apj;
287   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
288   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
289 
290   /* malloc apa to store dense row A[i,:]*P */
291   ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
292   ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
293 
294   ptap->apa = apa;
295 
296   /* create and assemble symbolic parallel matrix Cmpi */
297   /*----------------------------------------------------*/
298   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
299   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
300   ierr = MatSetBlockSizes(Cmpi,A->rmap->bs,P->cmap->bs);CHKERRQ(ierr);
301 
302   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
303   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
304   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
305   for (i=0; i<am; i++) {
306     row  = i + rstart;
307     apnz = api[i+1] - api[i];
308     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
309     apj += apnz;
310   }
311   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
312   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
313 
314   ptap->destroy        = Cmpi->ops->destroy;
315   ptap->duplicate      = Cmpi->ops->duplicate;
316   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
317   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
318 
319   /* attach the supporting struct to Cmpi for reuse */
320   c       = (Mat_MPIAIJ*)Cmpi->data;
321   c->ptap = ptap;
322 
323   *C = Cmpi;
324 
325   /* set MatInfo */
326   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
327   if (afill < 1.0) afill = 1.0;
328   Cmpi->info.mallocs           = nspacedouble;
329   Cmpi->info.fill_ratio_given  = fill;
330   Cmpi->info.fill_ratio_needed = afill;
331 
332 #if defined(PETSC_USE_INFO)
333   if (api[am]) {
334     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
335     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
336   } else {
337     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
338   }
339 #endif
340   PetscFunctionReturn(0);
341 }
342 
343 #undef __FUNCT__
344 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
345 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
346 {
347   PetscErrorCode ierr;
348 
349   PetscFunctionBegin;
350   if (scall == MAT_INITIAL_MATRIX) {
351     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
352   }
353   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
354   PetscFunctionReturn(0);
355 }
356 
357 typedef struct {
358   Mat         workB;
359   PetscScalar *rvalues,*svalues;
360   MPI_Request *rwaits,*swaits;
361 } MPIAIJ_MPIDense;
362 
363 #undef __FUNCT__
364 #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
365 PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
366 {
367   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
368   PetscErrorCode  ierr;
369 
370   PetscFunctionBegin;
371   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
372   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
373   ierr = PetscFree(contents);CHKERRQ(ierr);
374   PetscFunctionReturn(0);
375 }
376 
377 #undef __FUNCT__
378 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
379 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
380 {
381   PetscErrorCode         ierr;
382   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
383   PetscInt               nz   = aij->B->cmap->n;
384   PetscContainer         container;
385   MPIAIJ_MPIDense        *contents;
386   VecScatter             ctx   = aij->Mvctx;
387   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
388   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
389   PetscInt               m     = A->rmap->n,n=B->cmap->n;
390 
391   PetscFunctionBegin;
392   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
393   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
394   ierr = MatSetBlockSizes(*C,A->rmap->bs,B->cmap->bs);CHKERRQ(ierr);
395   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
396   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
397   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
398   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
399 
400   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
401 
402   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
403   /* Create work matrix used to store off processor rows of B needed for local product */
404   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
405   /* Create work arrays needed */
406   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
407                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
408                       from->n,MPI_Request,&contents->rwaits,
409                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
410 
411   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
412   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
413   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
414   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
415   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
416   PetscFunctionReturn(0);
417 }
418 
419 #undef __FUNCT__
420 #define __FUNCT__ "MatMPIDenseScatter"
421 /*
422     Performs an efficient scatter on the rows of B needed by this process; this is
423     a modification of the VecScatterBegin_() routines.
424 */
425 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
426 {
427   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
428   PetscErrorCode         ierr;
429   PetscScalar            *b,*w,*svalues,*rvalues;
430   VecScatter             ctx   = aij->Mvctx;
431   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
432   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
433   PetscInt               i,j,k;
434   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
435   PetscMPIInt            *sprocs,*rprocs,nrecvs;
436   MPI_Request            *swaits,*rwaits;
437   MPI_Comm               comm;
438   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
439   MPI_Status             status;
440   MPIAIJ_MPIDense        *contents;
441   PetscContainer         container;
442   Mat                    workB;
443 
444   PetscFunctionBegin;
445   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
446   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
447   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
448   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
449 
450   workB = *outworkB = contents->workB;
451   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
452   sindices = to->indices;
453   sstarts  = to->starts;
454   sprocs   = to->procs;
455   swaits   = contents->swaits;
456   svalues  = contents->svalues;
457 
458   rindices = from->indices;
459   rstarts  = from->starts;
460   rprocs   = from->procs;
461   rwaits   = contents->rwaits;
462   rvalues  = contents->rvalues;
463 
464   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
465   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
466 
467   for (i=0; i<from->n; i++) {
468     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
469   }
470 
471   for (i=0; i<to->n; i++) {
472     /* pack a message at a time */
473     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
474       for (k=0; k<ncols; k++) {
475         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
476       }
477     }
478     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
479   }
480 
481   nrecvs = from->n;
482   while (nrecvs) {
483     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
484     nrecvs--;
485     /* unpack a message at a time */
486     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
487       for (k=0; k<ncols; k++) {
488         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
489       }
490     }
491   }
492   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
493 
494   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
495   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
496   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
497   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
501 
502 #undef __FUNCT__
503 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
504 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
505 {
506   PetscErrorCode ierr;
507   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
508   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
509   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
510   Mat            workB;
511 
512   PetscFunctionBegin;
513   /* diagonal block of A times all local rows of B*/
514   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
515 
516   /* get off processor parts of B needed to complete the product */
517   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
518 
519   /* off-diagonal block of A times nonlocal rows of B */
520   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
521   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
522   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
523   PetscFunctionReturn(0);
524 }
525 
526 #undef __FUNCT__
527 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
528 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C)
529 {
530   PetscErrorCode ierr;
531   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
532   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
533   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
534   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
535   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
536   Mat_SeqAIJ     *p_loc,*p_oth;
537   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
538   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
539   PetscInt       cm          = C->rmap->n,anz,pnz;
540   Mat_PtAPMPI    *ptap       = c->ptap;
541   PetscScalar    *apa_sparse = ptap->apa;
542   PetscInt       *api,*apj,*apJ,i,j,k,row;
543   PetscInt       cstart = C->cmap->rstart;
544   PetscInt       cdnz,conz,k0,k1,nextp;
545 
546   PetscFunctionBegin;
547   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
548   /*-----------------------------------------------------*/
549   /* update numerical values of P_oth and P_loc */
550   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
551   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
552 
553   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
554   /*----------------------------------------------------------*/
555   /* get data from symbolic products */
556   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
557   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
558   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
559   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
560 
561   api = ptap->api;
562   apj = ptap->apj;
563   for (i=0; i<cm; i++) {
564     apJ = apj + api[i];
565 
566     /* diagonal portion of A */
567     anz = adi[i+1] - adi[i];
568     adj = ad->j + adi[i];
569     ada = ad->a + adi[i];
570     for (j=0; j<anz; j++) {
571       row = adj[j];
572       pnz = pi_loc[row+1] - pi_loc[row];
573       pj  = pj_loc + pi_loc[row];
574       pa  = pa_loc + pi_loc[row];
575       /* perform sparse axpy */
576       valtmp = ada[j];
577       nextp  = 0;
578       for (k=0; nextp<pnz; k++) {
579         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
580           apa_sparse[k] += valtmp*pa[nextp++];
581         }
582       }
583       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
584     }
585 
586     /* off-diagonal portion of A */
587     anz = aoi[i+1] - aoi[i];
588     aoj = ao->j + aoi[i];
589     aoa = ao->a + aoi[i];
590     for (j=0; j<anz; j++) {
591       row = aoj[j];
592       pnz = pi_oth[row+1] - pi_oth[row];
593       pj  = pj_oth + pi_oth[row];
594       pa  = pa_oth + pi_oth[row];
595       /* perform sparse axpy */
596       valtmp = aoa[j];
597       nextp  = 0;
598       for (k=0; nextp<pnz; k++) {
599         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
600           apa_sparse[k] += valtmp*pa[nextp++];
601         }
602       }
603       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
604     }
605 
606     /* set values in C */
607     cdnz = cd->i[i+1] - cd->i[i];
608     conz = co->i[i+1] - co->i[i];
609 
610     /* 1st off-diagoanl part of C */
611     ca = coa + co->i[i];
612     k  = 0;
613     for (k0=0; k0<conz; k0++) {
614       if (apJ[k] >= cstart) break;
615       ca[k0]        = apa_sparse[k];
616       apa_sparse[k] = 0.0;
617       k++;
618     }
619 
620     /* diagonal part of C */
621     ca = cda + cd->i[i];
622     for (k1=0; k1<cdnz; k1++) {
623       ca[k1]        = apa_sparse[k];
624       apa_sparse[k] = 0.0;
625       k++;
626     }
627 
628     /* 2nd off-diagoanl part of C */
629     ca = coa + co->i[i];
630     for (; k0<conz; k0++) {
631       ca[k0]        = apa_sparse[k];
632       apa_sparse[k] = 0.0;
633       k++;
634     }
635   }
636   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
637   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
638   PetscFunctionReturn(0);
639 }
640 
641 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */
642 #undef __FUNCT__
643 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
644 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C)
645 {
646   PetscErrorCode     ierr;
647   MPI_Comm           comm;
648   Mat                Cmpi;
649   Mat_PtAPMPI        *ptap;
650   PetscFreeSpaceList free_space = NULL,current_space=NULL;
651   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
652   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
653   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
654   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
655   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
656   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
657   PetscInt           nlnk_max,armax,prmax;
658   PetscReal          afill;
659   PetscScalar        *apa;
660 
661   PetscFunctionBegin;
662   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
663   /* create struct Mat_PtAPMPI and attached it to C later */
664   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
665 
666   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
667   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
668 
669   /* get P_loc by taking all local rows of P */
670   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
671 
672   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
673   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
674   pi_loc = p_loc->i; pj_loc = p_loc->j;
675   pi_oth = p_oth->i; pj_oth = p_oth->j;
676 
677   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
678   /*-------------------------------------------------------------------*/
679   ierr      = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
680   ptap->api = api;
681   api[0]    = 0;
682 
683   /* create and initialize a linked list */
684   armax    = ad->rmax+ao->rmax;
685   prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
686   nlnk_max = armax*prmax;
687   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
688   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
689 
690   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
691   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
692 
693   current_space = free_space;
694 
695   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
696   for (i=0; i<am; i++) {
697     /* diagonal portion of A */
698     nzi = adi[i+1] - adi[i];
699     for (j=0; j<nzi; j++) {
700       row  = *adj++;
701       pnz  = pi_loc[row+1] - pi_loc[row];
702       Jptr = pj_loc + pi_loc[row];
703       /* add non-zero cols of P into the sorted linked list lnk */
704       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
705     }
706     /* off-diagonal portion of A */
707     nzi = aoi[i+1] - aoi[i];
708     for (j=0; j<nzi; j++) {
709       row  = *aoj++;
710       pnz  = pi_oth[row+1] - pi_oth[row];
711       Jptr = pj_oth + pi_oth[row];
712       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
713     }
714 
715     apnz     = *lnk;
716     api[i+1] = api[i] + apnz;
717     if (apnz > apnz_max) apnz_max = apnz;
718 
719     /* if free space is not available, double the total space in the list */
720     if (current_space->local_remaining<apnz) {
721       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
722       nspacedouble++;
723     }
724 
725     /* Copy data into free space, then initialize lnk */
726     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
727     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
728 
729     current_space->array           += apnz;
730     current_space->local_used      += apnz;
731     current_space->local_remaining -= apnz;
732   }
733 
734   /* Allocate space for apj, initialize apj, and */
735   /* destroy list of free space and other temporary array(s) */
736   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
737   apj  = ptap->apj;
738   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
739   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
740 
741   /* create and assemble symbolic parallel matrix Cmpi */
742   /*----------------------------------------------------*/
743   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
744   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
745   ierr = MatSetBlockSizes(Cmpi,A->rmap->bs,P->cmap->bs);CHKERRQ(ierr);
746   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
747   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
748   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
749 
750   /* malloc apa for assembly Cmpi */
751   ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
752   ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr);
753 
754   ptap->apa = apa;
755   for (i=0; i<am; i++) {
756     row  = i + rstart;
757     apnz = api[i+1] - api[i];
758     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
759     apj += apnz;
760   }
761   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
762   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
763 
764   ptap->destroy             = Cmpi->ops->destroy;
765   ptap->duplicate           = Cmpi->ops->duplicate;
766   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
767   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
768   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
769 
770   /* attach the supporting struct to Cmpi for reuse */
771   c       = (Mat_MPIAIJ*)Cmpi->data;
772   c->ptap = ptap;
773 
774   *C = Cmpi;
775 
776   /* set MatInfo */
777   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
778   if (afill < 1.0) afill = 1.0;
779   Cmpi->info.mallocs           = nspacedouble;
780   Cmpi->info.fill_ratio_given  = fill;
781   Cmpi->info.fill_ratio_needed = afill;
782 
783 #if defined(PETSC_USE_INFO)
784   if (api[am]) {
785     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
786     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
787   } else {
788     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
789   }
790 #endif
791   PetscFunctionReturn(0);
792 }
793 
794 /*-------------------------------------------------------------------------*/
795 #undef __FUNCT__
796 #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
797 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
798 {
799   PetscErrorCode ierr;
800   PetscBool      scalable=PETSC_TRUE,viamatmatmult=PETSC_FALSE;
801 
802   PetscFunctionBegin;
803   ierr = PetscOptionsBool("-mattransposematmult_viamatmatmult","Use R=Pt and C=R*A","",viamatmatmult,&viamatmatmult,NULL);CHKERRQ(ierr);
804   if (viamatmatmult) {
805     Mat         Pt;
806     Mat_PtAPMPI *ptap;
807     Mat_MPIAIJ  *c;
808     if (scall == MAT_INITIAL_MATRIX) {
809       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
810       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
811 
812       c        = (Mat_MPIAIJ*)(*C)->data;
813       ptap     = c->ptap;
814       ptap->Pt = Pt;
815     } else if (scall == MAT_REUSE_MATRIX) {
816       c    = (Mat_MPIAIJ*)(*C)->data;
817       ptap = c->ptap;
818       Pt   = ptap->Pt;
819       ierr = MatTranspose(P,scall,&Pt);CHKERRQ(ierr);
820       ierr = MatMatMult(Pt,A,scall,fill,C);CHKERRQ(ierr);
821     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not supported");
822     PetscFunctionReturn(0);
823   }
824 
825   if (scall == MAT_INITIAL_MATRIX) {
826     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
827     ierr = PetscOptionsBool("-mattransposematmult_scalable","Use a scalable but slower C=Pt*A","",scalable,&scalable,NULL);CHKERRQ(ierr);
828     if  (scalable) {
829       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(P,A,fill,C);CHKERRQ(ierr);
830     } else {
831       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
832     }
833     ierr = PetscOptionsEnd();CHKERRQ(ierr);
834   }
835   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
836   PetscFunctionReturn(0);
837 }
838 
839 #undef __FUNCT__
840 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
841 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
842 {
843   PetscErrorCode      ierr;
844   Mat_Merge_SeqsToMPI *merge;
845   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
846   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
847   Mat_PtAPMPI         *ptap;
848   PetscInt            *adj,*aJ;
849   PetscInt            i,j,k,anz,pnz,row,*cj;
850   MatScalar           *ada,*aval,*ca,valtmp;
851   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
852   MPI_Comm            comm;
853   PetscMPIInt         size,rank,taga,*len_s;
854   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
855   PetscInt            **buf_ri,**buf_rj;
856   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
857   MPI_Request         *s_waits,*r_waits;
858   MPI_Status          *status;
859   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
860   PetscInt            *ai,*aj,*coi,*coj;
861   PetscInt            *poJ,*pdJ;
862   Mat                 A_loc;
863   Mat_SeqAIJ          *a_loc;
864 
865   PetscFunctionBegin;
866   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
867   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
868   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
869 
870   ptap  = c->ptap;
871   merge = ptap->merge;
872 
873   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
874   /*--------------------------------------------------------------*/
875   /* get data from symbolic products */
876   coi  = merge->coi; coj = merge->coj;
877   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
878   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
879 
880   bi     = merge->bi; bj = merge->bj;
881   owners = merge->rowmap->range;
882   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
883   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
884 
885   /* get A_loc by taking all local rows of A */
886   A_loc = ptap->A_loc;
887   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
888   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
889   ai    = a_loc->i;
890   aj    = a_loc->j;
891 
892   ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&aval);CHKERRQ(ierr); /* non-scalable!!! */
893   ierr = PetscMemzero(aval,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
894 
895   for (i=0; i<am; i++) {
896     /* 2-a) put A[i,:] to dense array aval */
897     anz = ai[i+1] - ai[i];
898     adj = aj + ai[i];
899     ada = a_loc->a + ai[i];
900     for (j=0; j<anz; j++) {
901       aval[adj[j]] = ada[j];
902     }
903 
904     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
905     /*--------------------------------------------------------------*/
906     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
907     pnz = po->i[i+1] - po->i[i];
908     poJ = po->j + po->i[i];
909     pA  = po->a + po->i[i];
910     for (j=0; j<pnz; j++) {
911       row = poJ[j];
912       cnz = coi[row+1] - coi[row];
913       cj  = coj + coi[row];
914       ca  = coa + coi[row];
915       /* perform dense axpy */
916       valtmp = pA[j];
917       for (k=0; k<cnz; k++) {
918         ca[k] += valtmp*aval[cj[k]];
919       }
920       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
921     }
922 
923     /* put the value into Cd (diagonal part) */
924     pnz = pd->i[i+1] - pd->i[i];
925     pdJ = pd->j + pd->i[i];
926     pA  = pd->a + pd->i[i];
927     for (j=0; j<pnz; j++) {
928       row = pdJ[j];
929       cnz = bi[row+1] - bi[row];
930       cj  = bj + bi[row];
931       ca  = ba + bi[row];
932       /* perform dense axpy */
933       valtmp = pA[j];
934       for (k=0; k<cnz; k++) {
935         ca[k] += valtmp*aval[cj[k]];
936       }
937       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
938     }
939 
940     /* zero the current row of Pt*A */
941     aJ = aj + ai[i];
942     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
943   }
944 
945   /* 3) send and recv matrix values coa */
946   /*------------------------------------*/
947   buf_ri = merge->buf_ri;
948   buf_rj = merge->buf_rj;
949   len_s  = merge->len_s;
950   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
951   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
952 
953   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
954   for (proc=0,k=0; proc<size; proc++) {
955     if (!len_s[proc]) continue;
956     i    = merge->owners_co[proc];
957     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
958     k++;
959   }
960   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
961   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
962 
963   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
964   ierr = PetscFree(r_waits);CHKERRQ(ierr);
965   ierr = PetscFree(coa);CHKERRQ(ierr);
966 
967   /* 4) insert local Cseq and received values into Cmpi */
968   /*----------------------------------------------------*/
969   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
970   for (k=0; k<merge->nrecv; k++) {
971     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
972     nrows       = *(buf_ri_k[k]);
973     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
974     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
975   }
976 
977   for (i=0; i<cm; i++) {
978     row  = owners[rank] + i; /* global row index of C_seq */
979     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
980     ba_i = ba + bi[i];
981     bnz  = bi[i+1] - bi[i];
982     /* add received vals into ba */
983     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
984       /* i-th row */
985       if (i == *nextrow[k]) {
986         cnz    = *(nextci[k]+1) - *nextci[k];
987         cj     = buf_rj[k] + *(nextci[k]);
988         ca     = abuf_r[k] + *(nextci[k]);
989         nextcj = 0;
990         for (j=0; nextcj<cnz; j++) {
991           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
992             ba_i[j] += ca[nextcj++];
993           }
994         }
995         nextrow[k]++; nextci[k]++;
996         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
997       }
998     }
999     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1000   }
1001   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1002   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1003 
1004   ierr = PetscFree(ba);CHKERRQ(ierr);
1005   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1006   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1007   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1008   ierr = PetscFree(aval);CHKERRQ(ierr);
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1013 #undef __FUNCT__
1014 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
1015 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1016 {
1017   PetscErrorCode      ierr;
1018   Mat                 Cmpi,A_loc,POt,PDt;
1019   Mat_PtAPMPI         *ptap;
1020   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1021   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1022   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1023   PetscInt            nnz;
1024   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1025   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1026   PetscBT             lnkbt;
1027   MPI_Comm            comm;
1028   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1029   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1030   PetscInt            len,proc,*dnz,*onz,*owners;
1031   PetscInt            nzi,*bi,*bj;
1032   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1033   MPI_Request         *swaits,*rwaits;
1034   MPI_Status          *sstatus,rstatus;
1035   Mat_Merge_SeqsToMPI *merge;
1036   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1037   PetscReal           afill  =1.0,afill_tmp;
1038   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1039   PetscScalar         *vals;
1040   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1041 
1042   PetscFunctionBegin;
1043   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1044   /* check if matrix local sizes are compatible */
1045   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1046     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1047   }
1048 
1049   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1050   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1051 
1052   /* create struct Mat_PtAPMPI and attached it to C later */
1053   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1054 
1055   /* get A_loc by taking all local rows of A */
1056   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1057 
1058   ptap->A_loc = A_loc;
1059 
1060   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1061   ai    = a_loc->i;
1062   aj    = a_loc->j;
1063 
1064   /* determine symbolic Co=(p->B)^T*A - send to others */
1065   /*----------------------------------------------------*/
1066   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1067   pdt  = (Mat_SeqAIJ*)PDt->data;
1068   pdti = pdt->i; pdtj = pdt->j;
1069 
1070   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1071   pot  = (Mat_SeqAIJ*)POt->data;
1072   poti = pot->i; potj = pot->j;
1073 
1074   /* then, compute symbolic Co = (p->B)^T*A */
1075   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1076   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1077   coi[0] = 0;
1078 
1079   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1080   nnz           = fill*(poti[pon] + ai[am]);
1081   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1082   current_space = free_space;
1083 
1084   /* create and initialize a linked list */
1085   i     = PetscMax(pdt->rmax,pot->rmax);
1086   Crmax = i*a_loc->rmax*size;
1087   if (!Crmax || Crmax > aN) Crmax = aN;
1088   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1089 
1090   for (i=0; i<pon; i++) {
1091     pnz = poti[i+1] - poti[i];
1092     ptJ = potj + poti[i];
1093     for (j=0; j<pnz; j++) {
1094       row  = ptJ[j]; /* row of A_loc == col of Pot */
1095       anz  = ai[row+1] - ai[row];
1096       Jptr = aj + ai[row];
1097       /* add non-zero cols of AP into the sorted linked list lnk */
1098       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1099     }
1100     nnz = lnk[0];
1101 
1102     /* If free space is not available, double the total space in the list */
1103     if (current_space->local_remaining<nnz) {
1104       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1105       nspacedouble++;
1106     }
1107 
1108     /* Copy data into free space, and zero out denserows */
1109     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1110 
1111     current_space->array           += nnz;
1112     current_space->local_used      += nnz;
1113     current_space->local_remaining -= nnz;
1114 
1115     coi[i+1] = coi[i] + nnz;
1116   }
1117 
1118   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1119   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1120 
1121   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1122   if (afill_tmp > afill) afill = afill_tmp;
1123 
1124   /* send j-array (coj) of Co to other processors */
1125   /*----------------------------------------------*/
1126   /* determine row ownership */
1127   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1128   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1129 
1130   merge->rowmap->n  = pn;
1131   merge->rowmap->bs = 1;
1132 
1133   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1134   owners = merge->rowmap->range;
1135 
1136   /* determine the number of messages to send, their lengths */
1137   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1138   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1139   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1140 
1141   len_s        = merge->len_s;
1142   merge->nsend = 0;
1143 
1144   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1145   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1146 
1147   proc = 0;
1148   for (i=0; i<pon; i++) {
1149     while (prmap[i] >= owners[proc+1]) proc++;
1150     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1151     len_s[proc] += coi[i+1] - coi[i];
1152   }
1153 
1154   len          = 0; /* max length of buf_si[] */
1155   owners_co[0] = 0;
1156   for (proc=0; proc<size; proc++) {
1157     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1158     if (len_si[proc]) {
1159       merge->nsend++;
1160       len_si[proc] = 2*(len_si[proc] + 1);
1161       len         += len_si[proc];
1162     }
1163   }
1164 
1165   /* determine the number and length of messages to receive for coi and coj  */
1166   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1167   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1168 
1169   /* post the Irecv and Isend of coj */
1170   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1171   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1172   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1173   for (proc=0, k=0; proc<size; proc++) {
1174     if (!len_s[proc]) continue;
1175     i    = owners_co[proc];
1176     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1177     k++;
1178   }
1179 
1180   /* receives and sends of coj are complete */
1181   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1182   for (i=0; i<merge->nrecv; i++) {
1183     PetscMPIInt icompleted;
1184     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1185   }
1186   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1187   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1188 
1189   /* send and recv coi */
1190   /*-------------------*/
1191   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1192   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1193   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1194   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1195   for (proc=0,k=0; proc<size; proc++) {
1196     if (!len_s[proc]) continue;
1197     /* form outgoing message for i-structure:
1198          buf_si[0]:                 nrows to be sent
1199                [1:nrows]:           row index (global)
1200                [nrows+1:2*nrows+1]: i-structure index
1201     */
1202     /*-------------------------------------------*/
1203     nrows       = len_si[proc]/2 - 1;
1204     buf_si_i    = buf_si + nrows+1;
1205     buf_si[0]   = nrows;
1206     buf_si_i[0] = 0;
1207     nrows       = 0;
1208     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1209       nzi               = coi[i+1] - coi[i];
1210       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1211       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1212       nrows++;
1213     }
1214     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1215     k++;
1216     buf_si += len_si[proc];
1217   }
1218   i = merge->nrecv;
1219   while (i--) {
1220     PetscMPIInt icompleted;
1221     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1222   }
1223   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1224   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1225   ierr = PetscFree(len_si);CHKERRQ(ierr);
1226   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1227   ierr = PetscFree(swaits);CHKERRQ(ierr);
1228   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1229   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1230 
1231   /* compute the local portion of C (mpi mat) */
1232   /*------------------------------------------*/
1233   /* allocate bi array and free space for accumulating nonzero column info */
1234   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1235   bi[0] = 0;
1236 
1237   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1238   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1239   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1240   current_space = free_space;
1241 
1242   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1243   for (k=0; k<merge->nrecv; k++) {
1244     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1245     nrows       = *buf_ri_k[k];
1246     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1247     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1248   }
1249 
1250   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1251   rmax = 0;
1252   for (i=0; i<pn; i++) {
1253     /* add pdt[i,:]*AP into lnk */
1254     pnz = pdti[i+1] - pdti[i];
1255     ptJ = pdtj + pdti[i];
1256     for (j=0; j<pnz; j++) {
1257       row  = ptJ[j];  /* row of AP == col of Pt */
1258       anz  = ai[row+1] - ai[row];
1259       Jptr = aj + ai[row];
1260       /* add non-zero cols of AP into the sorted linked list lnk */
1261       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1262     }
1263 
1264     /* add received col data into lnk */
1265     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1266       if (i == *nextrow[k]) { /* i-th row */
1267         nzi  = *(nextci[k]+1) - *nextci[k];
1268         Jptr = buf_rj[k] + *nextci[k];
1269         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1270         nextrow[k]++; nextci[k]++;
1271       }
1272     }
1273     nnz = lnk[0];
1274 
1275     /* if free space is not available, make more free space */
1276     if (current_space->local_remaining<nnz) {
1277       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1278       nspacedouble++;
1279     }
1280     /* copy data into free space, then initialize lnk */
1281     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1282     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1283 
1284     current_space->array           += nnz;
1285     current_space->local_used      += nnz;
1286     current_space->local_remaining -= nnz;
1287 
1288     bi[i+1] = bi[i] + nnz;
1289     if (nnz > rmax) rmax = nnz;
1290   }
1291   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1292 
1293   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1294   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1295 
1296   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1297   if (afill_tmp > afill) afill = afill_tmp;
1298   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
1299   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1300   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1301 
1302   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1303   /*----------------------------------------------------------------------------------*/
1304   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1305   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1306 
1307   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1308   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1309   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1310   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1311   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1312   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1313   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1314   for (i=0; i<pn; i++) {
1315     row  = i + rstart;
1316     nnz  = bi[i+1] - bi[i];
1317     Jptr = bj + bi[i];
1318     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1319   }
1320   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1321   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1322   ierr = PetscFree(vals);CHKERRQ(ierr);
1323 
1324   merge->bi        = bi;
1325   merge->bj        = bj;
1326   merge->coi       = coi;
1327   merge->coj       = coj;
1328   merge->buf_ri    = buf_ri;
1329   merge->buf_rj    = buf_rj;
1330   merge->owners_co = owners_co;
1331   merge->destroy   = Cmpi->ops->destroy;
1332   merge->duplicate = Cmpi->ops->duplicate;
1333 
1334   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1335   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1336 
1337   /* attach the supporting struct to Cmpi for reuse */
1338   c           = (Mat_MPIAIJ*)Cmpi->data;
1339   c->ptap     = ptap;
1340   ptap->api   = NULL;
1341   ptap->apj   = NULL;
1342   ptap->merge = merge;
1343   ptap->rmax  = rmax;
1344 
1345   *C = Cmpi;
1346 #if defined(PETSC_USE_INFO)
1347   if (bi[pn] != 0) {
1348     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1349     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1350   } else {
1351     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1352   }
1353 #endif
1354   PetscFunctionReturn(0);
1355 }
1356 
1357 #undef __FUNCT__
1358 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
1359 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,Mat C)
1360 {
1361   PetscErrorCode      ierr;
1362   Mat_Merge_SeqsToMPI *merge;
1363   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1364   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1365   Mat_PtAPMPI         *ptap;
1366   PetscInt            *adj;
1367   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1368   MatScalar           *ada,*ca,valtmp;
1369   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1370   MPI_Comm            comm;
1371   PetscMPIInt         size,rank,taga,*len_s;
1372   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1373   PetscInt            **buf_ri,**buf_rj;
1374   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1375   MPI_Request         *s_waits,*r_waits;
1376   MPI_Status          *status;
1377   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1378   PetscInt            *ai,*aj,*coi,*coj;
1379   PetscInt            *poJ,*pdJ;
1380   Mat                 A_loc;
1381   Mat_SeqAIJ          *a_loc;
1382 
1383   PetscFunctionBegin;
1384   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1385   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1386   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1387 
1388   ptap  = c->ptap;
1389   merge = ptap->merge;
1390 
1391   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1392   /*------------------------------------------*/
1393   /* get data from symbolic products */
1394   coi    = merge->coi; coj = merge->coj;
1395   ierr   = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
1396   ierr   = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
1397   bi     = merge->bi; bj = merge->bj;
1398   owners = merge->rowmap->range;
1399   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1400   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1401 
1402   /* get A_loc by taking all local rows of A */
1403   A_loc = ptap->A_loc;
1404   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1405   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1406   ai    = a_loc->i;
1407   aj    = a_loc->j;
1408 
1409   for (i=0; i<am; i++) {
1410     anz = ai[i+1] - ai[i];
1411     adj = aj + ai[i];
1412     ada = a_loc->a + ai[i];
1413 
1414     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1415     /*-------------------------------------------------------------*/
1416     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1417     pnz = po->i[i+1] - po->i[i];
1418     poJ = po->j + po->i[i];
1419     pA  = po->a + po->i[i];
1420     for (j=0; j<pnz; j++) {
1421       row = poJ[j];
1422       cj  = coj + coi[row];
1423       ca  = coa + coi[row];
1424       /* perform sparse axpy */
1425       nexta  = 0;
1426       valtmp = pA[j];
1427       for (k=0; nexta<anz; k++) {
1428         if (cj[k] == adj[nexta]) {
1429           ca[k] += valtmp*ada[nexta];
1430           nexta++;
1431         }
1432       }
1433       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1434     }
1435 
1436     /* put the value into Cd (diagonal part) */
1437     pnz = pd->i[i+1] - pd->i[i];
1438     pdJ = pd->j + pd->i[i];
1439     pA  = pd->a + pd->i[i];
1440     for (j=0; j<pnz; j++) {
1441       row = pdJ[j];
1442       cj  = bj + bi[row];
1443       ca  = ba + bi[row];
1444       /* perform sparse axpy */
1445       nexta  = 0;
1446       valtmp = pA[j];
1447       for (k=0; nexta<anz; k++) {
1448         if (cj[k] == adj[nexta]) {
1449           ca[k] += valtmp*ada[nexta];
1450           nexta++;
1451         }
1452       }
1453       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1454     }
1455   }
1456 
1457   /* 3) send and recv matrix values coa */
1458   /*------------------------------------*/
1459   buf_ri = merge->buf_ri;
1460   buf_rj = merge->buf_rj;
1461   len_s  = merge->len_s;
1462   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1463   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1464 
1465   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
1466   for (proc=0,k=0; proc<size; proc++) {
1467     if (!len_s[proc]) continue;
1468     i    = merge->owners_co[proc];
1469     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1470     k++;
1471   }
1472   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1473   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1474 
1475   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1476   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1477   ierr = PetscFree(coa);CHKERRQ(ierr);
1478 
1479   /* 4) insert local Cseq and received values into Cmpi */
1480   /*----------------------------------------------------*/
1481   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1482   for (k=0; k<merge->nrecv; k++) {
1483     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1484     nrows       = *(buf_ri_k[k]);
1485     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1486     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1487   }
1488 
1489   for (i=0; i<cm; i++) {
1490     row  = owners[rank] + i; /* global row index of C_seq */
1491     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1492     ba_i = ba + bi[i];
1493     bnz  = bi[i+1] - bi[i];
1494     /* add received vals into ba */
1495     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1496       /* i-th row */
1497       if (i == *nextrow[k]) {
1498         cnz    = *(nextci[k]+1) - *nextci[k];
1499         cj     = buf_rj[k] + *(nextci[k]);
1500         ca     = abuf_r[k] + *(nextci[k]);
1501         nextcj = 0;
1502         for (j=0; nextcj<cnz; j++) {
1503           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1504             ba_i[j] += ca[nextcj++];
1505           }
1506         }
1507         nextrow[k]++; nextci[k]++;
1508         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1509       }
1510     }
1511     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1512   }
1513   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1514   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1515 
1516   ierr = PetscFree(ba);CHKERRQ(ierr);
1517   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1518   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1519   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1520   PetscFunctionReturn(0);
1521 }
1522 
1523 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1524 #undef __FUNCT__
1525 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
1526 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,PetscReal fill,Mat *C)
1527 {
1528   PetscErrorCode      ierr;
1529   Mat                 Cmpi,A_loc,POt,PDt;
1530   Mat_PtAPMPI         *ptap;
1531   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1532   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1533   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1534   PetscInt            nnz;
1535   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1536   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1537   MPI_Comm            comm;
1538   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1539   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1540   PetscInt            len,proc,*dnz,*onz,*owners;
1541   PetscInt            nzi,*bi,*bj;
1542   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1543   MPI_Request         *swaits,*rwaits;
1544   MPI_Status          *sstatus,rstatus;
1545   Mat_Merge_SeqsToMPI *merge;
1546   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1547   PetscReal           afill  =1.0,afill_tmp;
1548   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1549   PetscScalar         *vals;
1550   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1551 
1552   PetscFunctionBegin;
1553   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1554   /* check if matrix local sizes are compatible */
1555   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1556     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1557   }
1558 
1559   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1560   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1561 
1562   /* create struct Mat_PtAPMPI and attached it to C later */
1563   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1564 
1565   /* get A_loc by taking all local rows of A */
1566   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1567 
1568   ptap->A_loc = A_loc;
1569   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1570   ai          = a_loc->i;
1571   aj          = a_loc->j;
1572 
1573   /* determine symbolic Co=(p->B)^T*A - send to others */
1574   /*----------------------------------------------------*/
1575   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1576   pdt  = (Mat_SeqAIJ*)PDt->data;
1577   pdti = pdt->i; pdtj = pdt->j;
1578 
1579   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1580   pot  = (Mat_SeqAIJ*)POt->data;
1581   poti = pot->i; potj = pot->j;
1582 
1583   /* then, compute symbolic Co = (p->B)^T*A */
1584   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1585                          >= (num of nonzero rows of C_seq) - pn */
1586   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1587   coi[0] = 0;
1588 
1589   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1590   nnz           = fill*(poti[pon] + ai[am]);
1591   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1592   current_space = free_space;
1593 
1594   /* create and initialize a linked list */
1595   i     = PetscMax(pdt->rmax,pot->rmax);
1596   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1597   if (!Crmax || Crmax > aN) Crmax = aN;
1598   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1599 
1600   for (i=0; i<pon; i++) {
1601     pnz = poti[i+1] - poti[i];
1602     ptJ = potj + poti[i];
1603     for (j=0; j<pnz; j++) {
1604       row  = ptJ[j]; /* row of A_loc == col of Pot */
1605       anz  = ai[row+1] - ai[row];
1606       Jptr = aj + ai[row];
1607       /* add non-zero cols of AP into the sorted linked list lnk */
1608       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1609     }
1610     nnz = lnk[0];
1611 
1612     /* If free space is not available, double the total space in the list */
1613     if (current_space->local_remaining<nnz) {
1614       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1615       nspacedouble++;
1616     }
1617 
1618     /* Copy data into free space, and zero out denserows */
1619     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1620 
1621     current_space->array           += nnz;
1622     current_space->local_used      += nnz;
1623     current_space->local_remaining -= nnz;
1624 
1625     coi[i+1] = coi[i] + nnz;
1626   }
1627 
1628   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1629   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1630 
1631   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1632   if (afill_tmp > afill) afill = afill_tmp;
1633 
1634   /* send j-array (coj) of Co to other processors */
1635   /*----------------------------------------------*/
1636   /* determine row ownership */
1637   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1638   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1639 
1640   merge->rowmap->n  = pn;
1641   merge->rowmap->bs = 1;
1642 
1643   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1644   owners = merge->rowmap->range;
1645 
1646   /* determine the number of messages to send, their lengths */
1647   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1648   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1649   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1650 
1651   len_s        = merge->len_s;
1652   merge->nsend = 0;
1653 
1654   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1655   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1656 
1657   proc = 0;
1658   for (i=0; i<pon; i++) {
1659     while (prmap[i] >= owners[proc+1]) proc++;
1660     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1661     len_s[proc] += coi[i+1] - coi[i];
1662   }
1663 
1664   len          = 0; /* max length of buf_si[] */
1665   owners_co[0] = 0;
1666   for (proc=0; proc<size; proc++) {
1667     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1668     if (len_si[proc]) {
1669       merge->nsend++;
1670       len_si[proc] = 2*(len_si[proc] + 1);
1671       len         += len_si[proc];
1672     }
1673   }
1674 
1675   /* determine the number and length of messages to receive for coi and coj  */
1676   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1677   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1678 
1679   /* post the Irecv and Isend of coj */
1680   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1681   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1682   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1683   for (proc=0, k=0; proc<size; proc++) {
1684     if (!len_s[proc]) continue;
1685     i    = owners_co[proc];
1686     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1687     k++;
1688   }
1689 
1690   /* receives and sends of coj are complete */
1691   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1692   for (i=0; i<merge->nrecv; i++) {
1693     PetscMPIInt icompleted;
1694     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1695   }
1696   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1697   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1698 
1699   /* send and recv coi */
1700   /*-------------------*/
1701   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1702   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1703   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1704   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1705   for (proc=0,k=0; proc<size; proc++) {
1706     if (!len_s[proc]) continue;
1707     /* form outgoing message for i-structure:
1708          buf_si[0]:                 nrows to be sent
1709                [1:nrows]:           row index (global)
1710                [nrows+1:2*nrows+1]: i-structure index
1711     */
1712     /*-------------------------------------------*/
1713     nrows       = len_si[proc]/2 - 1;
1714     buf_si_i    = buf_si + nrows+1;
1715     buf_si[0]   = nrows;
1716     buf_si_i[0] = 0;
1717     nrows       = 0;
1718     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1719       nzi               = coi[i+1] - coi[i];
1720       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1721       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1722       nrows++;
1723     }
1724     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1725     k++;
1726     buf_si += len_si[proc];
1727   }
1728   i = merge->nrecv;
1729   while (i--) {
1730     PetscMPIInt icompleted;
1731     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1732   }
1733   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1734   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1735   ierr = PetscFree(len_si);CHKERRQ(ierr);
1736   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1737   ierr = PetscFree(swaits);CHKERRQ(ierr);
1738   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1739   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1740 
1741   /* compute the local portion of C (mpi mat) */
1742   /*------------------------------------------*/
1743   /* allocate bi array and free space for accumulating nonzero column info */
1744   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1745   bi[0] = 0;
1746 
1747   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1748   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1749   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1750   current_space = free_space;
1751 
1752   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1753   for (k=0; k<merge->nrecv; k++) {
1754     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1755     nrows       = *buf_ri_k[k];
1756     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1757     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1758   }
1759 
1760   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1761   rmax = 0;
1762   for (i=0; i<pn; i++) {
1763     /* add pdt[i,:]*AP into lnk */
1764     pnz = pdti[i+1] - pdti[i];
1765     ptJ = pdtj + pdti[i];
1766     for (j=0; j<pnz; j++) {
1767       row  = ptJ[j];  /* row of AP == col of Pt */
1768       anz  = ai[row+1] - ai[row];
1769       Jptr = aj + ai[row];
1770       /* add non-zero cols of AP into the sorted linked list lnk */
1771       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1772     }
1773 
1774     /* add received col data into lnk */
1775     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1776       if (i == *nextrow[k]) { /* i-th row */
1777         nzi  = *(nextci[k]+1) - *nextci[k];
1778         Jptr = buf_rj[k] + *nextci[k];
1779         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1780         nextrow[k]++; nextci[k]++;
1781       }
1782     }
1783     nnz = lnk[0];
1784 
1785     /* if free space is not available, make more free space */
1786     if (current_space->local_remaining<nnz) {
1787       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1788       nspacedouble++;
1789     }
1790     /* copy data into free space, then initialize lnk */
1791     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1792     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1793 
1794     current_space->array           += nnz;
1795     current_space->local_used      += nnz;
1796     current_space->local_remaining -= nnz;
1797 
1798     bi[i+1] = bi[i] + nnz;
1799     if (nnz > rmax) rmax = nnz;
1800   }
1801   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1802 
1803   ierr      = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1804   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1805   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1806   if (afill_tmp > afill) afill = afill_tmp;
1807   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1808   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1809   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1810 
1811   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1812   /*----------------------------------------------------------------------------------*/
1813   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1814   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1815 
1816   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1817   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1818   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1819   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1820   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1821   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1822   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1823   for (i=0; i<pn; i++) {
1824     row  = i + rstart;
1825     nnz  = bi[i+1] - bi[i];
1826     Jptr = bj + bi[i];
1827     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1828   }
1829   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1830   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1831   ierr = PetscFree(vals);CHKERRQ(ierr);
1832 
1833   merge->bi        = bi;
1834   merge->bj        = bj;
1835   merge->coi       = coi;
1836   merge->coj       = coj;
1837   merge->buf_ri    = buf_ri;
1838   merge->buf_rj    = buf_rj;
1839   merge->owners_co = owners_co;
1840   merge->destroy   = Cmpi->ops->destroy;
1841   merge->duplicate = Cmpi->ops->duplicate;
1842 
1843   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
1844   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1845 
1846   /* attach the supporting struct to Cmpi for reuse */
1847   c = (Mat_MPIAIJ*)Cmpi->data;
1848 
1849   c->ptap     = ptap;
1850   ptap->api   = NULL;
1851   ptap->apj   = NULL;
1852   ptap->merge = merge;
1853   ptap->rmax  = rmax;
1854   ptap->apa   = NULL;
1855 
1856   *C = Cmpi;
1857 #if defined(PETSC_USE_INFO)
1858   if (bi[pn] != 0) {
1859     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1860     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1861   } else {
1862     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1863   }
1864 #endif
1865   PetscFunctionReturn(0);
1866 }
1867