xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
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     CHKMEMQ;
474     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
475       for (k=0; k<ncols; k++) {
476         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
477       }
478     }
479     CHKMEMQ;
480     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
481   }
482 
483   nrecvs = from->n;
484   while (nrecvs) {
485     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
486     nrecvs--;
487     /* unpack a message at a time */
488     CHKMEMQ;
489     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
490       for (k=0; k<ncols; k++) {
491         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
492       }
493     }
494     CHKMEMQ;
495   }
496   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
497 
498   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
499   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
500   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
501   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
502   PetscFunctionReturn(0);
503 }
504 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
508 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
509 {
510   PetscErrorCode ierr;
511   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
512   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
513   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
514   Mat            workB;
515 
516   PetscFunctionBegin;
517   /* diagonal block of A times all local rows of B*/
518   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
519 
520   /* get off processor parts of B needed to complete the product */
521   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
522 
523   /* off-diagonal block of A times nonlocal rows of B */
524   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
525   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
526   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
527   PetscFunctionReturn(0);
528 }
529 
530 #undef __FUNCT__
531 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
532 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C)
533 {
534   PetscErrorCode ierr;
535   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
536   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
537   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
538   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
539   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
540   Mat_SeqAIJ     *p_loc,*p_oth;
541   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
542   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
543   PetscInt       cm          = C->rmap->n,anz,pnz;
544   Mat_PtAPMPI    *ptap       = c->ptap;
545   PetscScalar    *apa_sparse = ptap->apa;
546   PetscInt       *api,*apj,*apJ,i,j,k,row;
547   PetscInt       cstart = C->cmap->rstart;
548   PetscInt       cdnz,conz,k0,k1,nextp;
549 
550   PetscFunctionBegin;
551   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
552   /*-----------------------------------------------------*/
553   /* update numerical values of P_oth and P_loc */
554   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
555   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
556 
557   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
558   /*----------------------------------------------------------*/
559   /* get data from symbolic products */
560   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
561   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
562   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
563   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
564 
565   api = ptap->api;
566   apj = ptap->apj;
567   for (i=0; i<cm; i++) {
568     apJ = apj + api[i];
569 
570     /* diagonal portion of A */
571     anz = adi[i+1] - adi[i];
572     adj = ad->j + adi[i];
573     ada = ad->a + adi[i];
574     for (j=0; j<anz; j++) {
575       row = adj[j];
576       pnz = pi_loc[row+1] - pi_loc[row];
577       pj  = pj_loc + pi_loc[row];
578       pa  = pa_loc + pi_loc[row];
579       /* perform sparse axpy */
580       valtmp = ada[j];
581       nextp  = 0;
582       for (k=0; nextp<pnz; k++) {
583         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
584           apa_sparse[k] += valtmp*pa[nextp++];
585         }
586       }
587       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
588     }
589 
590     /* off-diagonal portion of A */
591     anz = aoi[i+1] - aoi[i];
592     aoj = ao->j + aoi[i];
593     aoa = ao->a + aoi[i];
594     for (j=0; j<anz; j++) {
595       row = aoj[j];
596       pnz = pi_oth[row+1] - pi_oth[row];
597       pj  = pj_oth + pi_oth[row];
598       pa  = pa_oth + pi_oth[row];
599       /* perform sparse axpy */
600       valtmp = aoa[j];
601       nextp  = 0;
602       for (k=0; nextp<pnz; k++) {
603         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
604           apa_sparse[k] += valtmp*pa[nextp++];
605         }
606       }
607       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
608     }
609 
610     /* set values in C */
611     cdnz = cd->i[i+1] - cd->i[i];
612     conz = co->i[i+1] - co->i[i];
613 
614     /* 1st off-diagoanl part of C */
615     ca = coa + co->i[i];
616     k  = 0;
617     for (k0=0; k0<conz; k0++) {
618       if (apJ[k] >= cstart) break;
619       ca[k0]        = apa_sparse[k];
620       apa_sparse[k] = 0.0;
621       k++;
622     }
623 
624     /* diagonal part of C */
625     ca = cda + cd->i[i];
626     for (k1=0; k1<cdnz; k1++) {
627       ca[k1]        = apa_sparse[k];
628       apa_sparse[k] = 0.0;
629       k++;
630     }
631 
632     /* 2nd off-diagoanl part of C */
633     ca = coa + co->i[i];
634     for (; k0<conz; k0++) {
635       ca[k0]        = apa_sparse[k];
636       apa_sparse[k] = 0.0;
637       k++;
638     }
639   }
640   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
641   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
642   PetscFunctionReturn(0);
643 }
644 
645 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */
646 #undef __FUNCT__
647 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
648 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C)
649 {
650   PetscErrorCode     ierr;
651   MPI_Comm           comm;
652   Mat                Cmpi;
653   Mat_PtAPMPI        *ptap;
654   PetscFreeSpaceList free_space = NULL,current_space=NULL;
655   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
656   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
657   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
658   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
659   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
660   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
661   PetscInt           nlnk_max,armax,prmax;
662   PetscReal          afill;
663   PetscScalar        *apa;
664 
665   PetscFunctionBegin;
666   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
667   /* create struct Mat_PtAPMPI and attached it to C later */
668   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
669 
670   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
671   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
672 
673   /* get P_loc by taking all local rows of P */
674   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
675 
676   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
677   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
678   pi_loc = p_loc->i; pj_loc = p_loc->j;
679   pi_oth = p_oth->i; pj_oth = p_oth->j;
680 
681   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
682   /*-------------------------------------------------------------------*/
683   ierr      = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
684   ptap->api = api;
685   api[0]    = 0;
686 
687   /* create and initialize a linked list */
688   armax    = ad->rmax+ao->rmax;
689   prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
690   nlnk_max = armax*prmax;
691   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
692   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
693 
694   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
695   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
696 
697   current_space = free_space;
698 
699   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
700   for (i=0; i<am; i++) {
701     /* diagonal portion of A */
702     nzi = adi[i+1] - adi[i];
703     for (j=0; j<nzi; j++) {
704       row  = *adj++;
705       pnz  = pi_loc[row+1] - pi_loc[row];
706       Jptr = pj_loc + pi_loc[row];
707       /* add non-zero cols of P into the sorted linked list lnk */
708       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
709     }
710     /* off-diagonal portion of A */
711     nzi = aoi[i+1] - aoi[i];
712     for (j=0; j<nzi; j++) {
713       row  = *aoj++;
714       pnz  = pi_oth[row+1] - pi_oth[row];
715       Jptr = pj_oth + pi_oth[row];
716       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
717     }
718 
719     apnz     = *lnk;
720     api[i+1] = api[i] + apnz;
721     if (apnz > apnz_max) apnz_max = apnz;
722 
723     /* if free space is not available, double the total space in the list */
724     if (current_space->local_remaining<apnz) {
725       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
726       nspacedouble++;
727     }
728 
729     /* Copy data into free space, then initialize lnk */
730     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
731     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
732 
733     current_space->array           += apnz;
734     current_space->local_used      += apnz;
735     current_space->local_remaining -= apnz;
736   }
737 
738   /* Allocate space for apj, initialize apj, and */
739   /* destroy list of free space and other temporary array(s) */
740   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
741   apj  = ptap->apj;
742   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
743   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
744 
745   /* create and assemble symbolic parallel matrix Cmpi */
746   /*----------------------------------------------------*/
747   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
748   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
749   ierr = MatSetBlockSizes(Cmpi,A->rmap->bs,P->cmap->bs);CHKERRQ(ierr);
750   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
751   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
752   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
753 
754   /* malloc apa for assembly Cmpi */
755   ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
756   ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr);
757 
758   ptap->apa = apa;
759   for (i=0; i<am; i++) {
760     row  = i + rstart;
761     apnz = api[i+1] - api[i];
762     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
763     apj += apnz;
764   }
765   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
766   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
767 
768   ptap->destroy             = Cmpi->ops->destroy;
769   ptap->duplicate           = Cmpi->ops->duplicate;
770   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
771   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
772   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
773 
774   /* attach the supporting struct to Cmpi for reuse */
775   c       = (Mat_MPIAIJ*)Cmpi->data;
776   c->ptap = ptap;
777 
778   *C = Cmpi;
779 
780   /* set MatInfo */
781   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
782   if (afill < 1.0) afill = 1.0;
783   Cmpi->info.mallocs           = nspacedouble;
784   Cmpi->info.fill_ratio_given  = fill;
785   Cmpi->info.fill_ratio_needed = afill;
786 
787 #if defined(PETSC_USE_INFO)
788   if (api[am]) {
789     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
790     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
791   } else {
792     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
793   }
794 #endif
795   PetscFunctionReturn(0);
796 }
797 
798 /*-------------------------------------------------------------------------*/
799 #undef __FUNCT__
800 #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
801 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
802 {
803   PetscErrorCode ierr;
804   PetscBool      scalable=PETSC_TRUE,viamatmatmult=PETSC_FALSE;
805 
806   PetscFunctionBegin;
807   ierr = PetscOptionsBool("-mattransposematmult_viamatmatmult","Use R=Pt and C=R*A","",viamatmatmult,&viamatmatmult,NULL);CHKERRQ(ierr);
808   if (viamatmatmult) {
809     Mat         Pt;
810     Mat_PtAPMPI *ptap;
811     Mat_MPIAIJ  *c;
812     if (scall == MAT_INITIAL_MATRIX) {
813       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
814       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
815 
816       c        = (Mat_MPIAIJ*)(*C)->data;
817       ptap     = c->ptap;
818       ptap->Pt = Pt;
819     } else if (scall == MAT_REUSE_MATRIX) {
820       c    = (Mat_MPIAIJ*)(*C)->data;
821       ptap = c->ptap;
822       Pt   = ptap->Pt;
823       ierr = MatTranspose(P,scall,&Pt);CHKERRQ(ierr);
824       ierr = MatMatMult(Pt,A,scall,fill,C);CHKERRQ(ierr);
825     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not supported");
826     PetscFunctionReturn(0);
827   }
828 
829   if (scall == MAT_INITIAL_MATRIX) {
830     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
831     ierr = PetscOptionsBool("-mattransposematmult_scalable","Use a scalable but slower C=Pt*A","",scalable,&scalable,NULL);CHKERRQ(ierr);
832     if  (scalable) {
833       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(P,A,fill,C);CHKERRQ(ierr);
834     } else {
835       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
836     }
837     ierr = PetscOptionsEnd();CHKERRQ(ierr);
838   }
839   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
840   PetscFunctionReturn(0);
841 }
842 
843 #undef __FUNCT__
844 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
845 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
846 {
847   PetscErrorCode      ierr;
848   Mat_Merge_SeqsToMPI *merge;
849   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
850   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
851   Mat_PtAPMPI         *ptap;
852   PetscInt            *adj,*aJ;
853   PetscInt            i,j,k,anz,pnz,row,*cj;
854   MatScalar           *ada,*aval,*ca,valtmp;
855   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
856   MPI_Comm            comm;
857   PetscMPIInt         size,rank,taga,*len_s;
858   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
859   PetscInt            **buf_ri,**buf_rj;
860   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
861   MPI_Request         *s_waits,*r_waits;
862   MPI_Status          *status;
863   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
864   PetscInt            *ai,*aj,*coi,*coj;
865   PetscInt            *poJ,*pdJ;
866   Mat                 A_loc;
867   Mat_SeqAIJ          *a_loc;
868 
869   PetscFunctionBegin;
870   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
871   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
872   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
873 
874   ptap  = c->ptap;
875   merge = ptap->merge;
876 
877   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
878   /*--------------------------------------------------------------*/
879   /* get data from symbolic products */
880   coi  = merge->coi; coj = merge->coj;
881   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
882   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
883 
884   bi     = merge->bi; bj = merge->bj;
885   owners = merge->rowmap->range;
886   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
887   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
888 
889   /* get A_loc by taking all local rows of A */
890   A_loc = ptap->A_loc;
891   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
892   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
893   ai    = a_loc->i;
894   aj    = a_loc->j;
895 
896   ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&aval);CHKERRQ(ierr); /* non-scalable!!! */
897   ierr = PetscMemzero(aval,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
898 
899   for (i=0; i<am; i++) {
900     /* 2-a) put A[i,:] to dense array aval */
901     anz = ai[i+1] - ai[i];
902     adj = aj + ai[i];
903     ada = a_loc->a + ai[i];
904     for (j=0; j<anz; j++) {
905       aval[adj[j]] = ada[j];
906     }
907 
908     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
909     /*--------------------------------------------------------------*/
910     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
911     pnz = po->i[i+1] - po->i[i];
912     poJ = po->j + po->i[i];
913     pA  = po->a + po->i[i];
914     for (j=0; j<pnz; j++) {
915       row = poJ[j];
916       cnz = coi[row+1] - coi[row];
917       cj  = coj + coi[row];
918       ca  = coa + coi[row];
919       /* perform dense axpy */
920       valtmp = pA[j];
921       for (k=0; k<cnz; k++) {
922         ca[k] += valtmp*aval[cj[k]];
923       }
924       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
925     }
926 
927     /* put the value into Cd (diagonal part) */
928     pnz = pd->i[i+1] - pd->i[i];
929     pdJ = pd->j + pd->i[i];
930     pA  = pd->a + pd->i[i];
931     for (j=0; j<pnz; j++) {
932       row = pdJ[j];
933       cnz = bi[row+1] - bi[row];
934       cj  = bj + bi[row];
935       ca  = ba + bi[row];
936       /* perform dense axpy */
937       valtmp = pA[j];
938       for (k=0; k<cnz; k++) {
939         ca[k] += valtmp*aval[cj[k]];
940       }
941       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
942     }
943 
944     /* zero the current row of Pt*A */
945     aJ = aj + ai[i];
946     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
947   }
948 
949   /* 3) send and recv matrix values coa */
950   /*------------------------------------*/
951   buf_ri = merge->buf_ri;
952   buf_rj = merge->buf_rj;
953   len_s  = merge->len_s;
954   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
955   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
956 
957   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
958   for (proc=0,k=0; proc<size; proc++) {
959     if (!len_s[proc]) continue;
960     i    = merge->owners_co[proc];
961     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
962     k++;
963   }
964   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
965   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
966 
967   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
968   ierr = PetscFree(r_waits);CHKERRQ(ierr);
969   ierr = PetscFree(coa);CHKERRQ(ierr);
970 
971   /* 4) insert local Cseq and received values into Cmpi */
972   /*----------------------------------------------------*/
973   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
974   for (k=0; k<merge->nrecv; k++) {
975     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
976     nrows       = *(buf_ri_k[k]);
977     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
978     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
979   }
980 
981   for (i=0; i<cm; i++) {
982     row  = owners[rank] + i; /* global row index of C_seq */
983     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
984     ba_i = ba + bi[i];
985     bnz  = bi[i+1] - bi[i];
986     /* add received vals into ba */
987     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
988       /* i-th row */
989       if (i == *nextrow[k]) {
990         cnz    = *(nextci[k]+1) - *nextci[k];
991         cj     = buf_rj[k] + *(nextci[k]);
992         ca     = abuf_r[k] + *(nextci[k]);
993         nextcj = 0;
994         for (j=0; nextcj<cnz; j++) {
995           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
996             ba_i[j] += ca[nextcj++];
997           }
998         }
999         nextrow[k]++; nextci[k]++;
1000         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1001       }
1002     }
1003     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1004   }
1005   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1006   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1007 
1008   ierr = PetscFree(ba);CHKERRQ(ierr);
1009   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1010   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1011   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1012   ierr = PetscFree(aval);CHKERRQ(ierr);
1013   PetscFunctionReturn(0);
1014 }
1015 
1016 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1017 #undef __FUNCT__
1018 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
1019 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1020 {
1021   PetscErrorCode      ierr;
1022   Mat                 Cmpi,A_loc,POt,PDt;
1023   Mat_PtAPMPI         *ptap;
1024   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1025   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1026   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1027   PetscInt            nnz;
1028   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1029   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1030   PetscBT             lnkbt;
1031   MPI_Comm            comm;
1032   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1033   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1034   PetscInt            len,proc,*dnz,*onz,*owners;
1035   PetscInt            nzi,*bi,*bj;
1036   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1037   MPI_Request         *swaits,*rwaits;
1038   MPI_Status          *sstatus,rstatus;
1039   Mat_Merge_SeqsToMPI *merge;
1040   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1041   PetscReal           afill  =1.0,afill_tmp;
1042   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1043   PetscScalar         *vals;
1044   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1045 
1046   PetscFunctionBegin;
1047   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1048   /* check if matrix local sizes are compatible */
1049   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1050     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);
1051   }
1052 
1053   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1054   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1055 
1056   /* create struct Mat_PtAPMPI and attached it to C later */
1057   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1058 
1059   /* get A_loc by taking all local rows of A */
1060   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1061 
1062   ptap->A_loc = A_loc;
1063 
1064   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1065   ai    = a_loc->i;
1066   aj    = a_loc->j;
1067 
1068   /* determine symbolic Co=(p->B)^T*A - send to others */
1069   /*----------------------------------------------------*/
1070   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1071   pdt  = (Mat_SeqAIJ*)PDt->data;
1072   pdti = pdt->i; pdtj = pdt->j;
1073 
1074   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1075   pot  = (Mat_SeqAIJ*)POt->data;
1076   poti = pot->i; potj = pot->j;
1077 
1078   /* then, compute symbolic Co = (p->B)^T*A */
1079   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1080   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1081   coi[0] = 0;
1082 
1083   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1084   nnz           = fill*(poti[pon] + ai[am]);
1085   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1086   current_space = free_space;
1087 
1088   /* create and initialize a linked list */
1089   i     = PetscMax(pdt->rmax,pot->rmax);
1090   Crmax = i*a_loc->rmax*size;
1091   if (!Crmax || Crmax > aN) Crmax = aN;
1092   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1093 
1094   for (i=0; i<pon; i++) {
1095     pnz = poti[i+1] - poti[i];
1096     ptJ = potj + poti[i];
1097     for (j=0; j<pnz; j++) {
1098       row  = ptJ[j]; /* row of A_loc == col of Pot */
1099       anz  = ai[row+1] - ai[row];
1100       Jptr = aj + ai[row];
1101       /* add non-zero cols of AP into the sorted linked list lnk */
1102       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1103     }
1104     nnz = lnk[0];
1105 
1106     /* If free space is not available, double the total space in the list */
1107     if (current_space->local_remaining<nnz) {
1108       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1109       nspacedouble++;
1110     }
1111 
1112     /* Copy data into free space, and zero out denserows */
1113     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1114 
1115     current_space->array           += nnz;
1116     current_space->local_used      += nnz;
1117     current_space->local_remaining -= nnz;
1118 
1119     coi[i+1] = coi[i] + nnz;
1120   }
1121 
1122   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1123   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1124 
1125   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1126   if (afill_tmp > afill) afill = afill_tmp;
1127 
1128   /* send j-array (coj) of Co to other processors */
1129   /*----------------------------------------------*/
1130   /* determine row ownership */
1131   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1132   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1133 
1134   merge->rowmap->n  = pn;
1135   merge->rowmap->bs = 1;
1136 
1137   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1138   owners = merge->rowmap->range;
1139 
1140   /* determine the number of messages to send, their lengths */
1141   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1142   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1143   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1144 
1145   len_s        = merge->len_s;
1146   merge->nsend = 0;
1147 
1148   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1149   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1150 
1151   proc = 0;
1152   for (i=0; i<pon; i++) {
1153     while (prmap[i] >= owners[proc+1]) proc++;
1154     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1155     len_s[proc] += coi[i+1] - coi[i];
1156   }
1157 
1158   len          = 0; /* max length of buf_si[] */
1159   owners_co[0] = 0;
1160   for (proc=0; proc<size; proc++) {
1161     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1162     if (len_si[proc]) {
1163       merge->nsend++;
1164       len_si[proc] = 2*(len_si[proc] + 1);
1165       len         += len_si[proc];
1166     }
1167   }
1168 
1169   /* determine the number and length of messages to receive for coi and coj  */
1170   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1171   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1172 
1173   /* post the Irecv and Isend of coj */
1174   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1175   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1176   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1177   for (proc=0, k=0; proc<size; proc++) {
1178     if (!len_s[proc]) continue;
1179     i    = owners_co[proc];
1180     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1181     k++;
1182   }
1183 
1184   /* receives and sends of coj are complete */
1185   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1186   for (i=0; i<merge->nrecv; i++) {
1187     PetscMPIInt icompleted;
1188     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1189   }
1190   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1191   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1192 
1193   /* send and recv coi */
1194   /*-------------------*/
1195   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1196   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1197   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1198   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1199   for (proc=0,k=0; proc<size; proc++) {
1200     if (!len_s[proc]) continue;
1201     /* form outgoing message for i-structure:
1202          buf_si[0]:                 nrows to be sent
1203                [1:nrows]:           row index (global)
1204                [nrows+1:2*nrows+1]: i-structure index
1205     */
1206     /*-------------------------------------------*/
1207     nrows       = len_si[proc]/2 - 1;
1208     buf_si_i    = buf_si + nrows+1;
1209     buf_si[0]   = nrows;
1210     buf_si_i[0] = 0;
1211     nrows       = 0;
1212     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1213       nzi               = coi[i+1] - coi[i];
1214       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1215       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1216       nrows++;
1217     }
1218     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1219     k++;
1220     buf_si += len_si[proc];
1221   }
1222   i = merge->nrecv;
1223   while (i--) {
1224     PetscMPIInt icompleted;
1225     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1226   }
1227   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1228   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1229   ierr = PetscFree(len_si);CHKERRQ(ierr);
1230   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1231   ierr = PetscFree(swaits);CHKERRQ(ierr);
1232   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1233   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1234 
1235   /* compute the local portion of C (mpi mat) */
1236   /*------------------------------------------*/
1237   /* allocate bi array and free space for accumulating nonzero column info */
1238   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1239   bi[0] = 0;
1240 
1241   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1242   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1243   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1244   current_space = free_space;
1245 
1246   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1247   for (k=0; k<merge->nrecv; k++) {
1248     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1249     nrows       = *buf_ri_k[k];
1250     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1251     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1252   }
1253 
1254   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1255   rmax = 0;
1256   for (i=0; i<pn; i++) {
1257     /* add pdt[i,:]*AP into lnk */
1258     pnz = pdti[i+1] - pdti[i];
1259     ptJ = pdtj + pdti[i];
1260     for (j=0; j<pnz; j++) {
1261       row  = ptJ[j];  /* row of AP == col of Pt */
1262       anz  = ai[row+1] - ai[row];
1263       Jptr = aj + ai[row];
1264       /* add non-zero cols of AP into the sorted linked list lnk */
1265       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1266     }
1267 
1268     /* add received col data into lnk */
1269     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1270       if (i == *nextrow[k]) { /* i-th row */
1271         nzi  = *(nextci[k]+1) - *nextci[k];
1272         Jptr = buf_rj[k] + *nextci[k];
1273         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1274         nextrow[k]++; nextci[k]++;
1275       }
1276     }
1277     nnz = lnk[0];
1278 
1279     /* if free space is not available, make more free space */
1280     if (current_space->local_remaining<nnz) {
1281       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1282       nspacedouble++;
1283     }
1284     /* copy data into free space, then initialize lnk */
1285     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1286     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1287 
1288     current_space->array           += nnz;
1289     current_space->local_used      += nnz;
1290     current_space->local_remaining -= nnz;
1291 
1292     bi[i+1] = bi[i] + nnz;
1293     if (nnz > rmax) rmax = nnz;
1294   }
1295   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1296 
1297   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1298   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1299 
1300   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1301   if (afill_tmp > afill) afill = afill_tmp;
1302   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
1303   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1304   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1305 
1306   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1307   /*----------------------------------------------------------------------------------*/
1308   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1309   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1310 
1311   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1312   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1313   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1314   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1315   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1316   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1317   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1318   for (i=0; i<pn; i++) {
1319     row  = i + rstart;
1320     nnz  = bi[i+1] - bi[i];
1321     Jptr = bj + bi[i];
1322     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1323   }
1324   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1325   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1326   ierr = PetscFree(vals);CHKERRQ(ierr);
1327 
1328   merge->bi        = bi;
1329   merge->bj        = bj;
1330   merge->coi       = coi;
1331   merge->coj       = coj;
1332   merge->buf_ri    = buf_ri;
1333   merge->buf_rj    = buf_rj;
1334   merge->owners_co = owners_co;
1335   merge->destroy   = Cmpi->ops->destroy;
1336   merge->duplicate = Cmpi->ops->duplicate;
1337 
1338   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1339   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1340 
1341   /* attach the supporting struct to Cmpi for reuse */
1342   c           = (Mat_MPIAIJ*)Cmpi->data;
1343   c->ptap     = ptap;
1344   ptap->api   = NULL;
1345   ptap->apj   = NULL;
1346   ptap->merge = merge;
1347   ptap->rmax  = rmax;
1348 
1349   *C = Cmpi;
1350 #if defined(PETSC_USE_INFO)
1351   if (bi[pn] != 0) {
1352     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1353     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1354   } else {
1355     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1356   }
1357 #endif
1358   PetscFunctionReturn(0);
1359 }
1360 
1361 #undef __FUNCT__
1362 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
1363 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,Mat C)
1364 {
1365   PetscErrorCode      ierr;
1366   Mat_Merge_SeqsToMPI *merge;
1367   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1368   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1369   Mat_PtAPMPI         *ptap;
1370   PetscInt            *adj;
1371   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1372   MatScalar           *ada,*ca,valtmp;
1373   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1374   MPI_Comm            comm;
1375   PetscMPIInt         size,rank,taga,*len_s;
1376   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1377   PetscInt            **buf_ri,**buf_rj;
1378   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1379   MPI_Request         *s_waits,*r_waits;
1380   MPI_Status          *status;
1381   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1382   PetscInt            *ai,*aj,*coi,*coj;
1383   PetscInt            *poJ,*pdJ;
1384   Mat                 A_loc;
1385   Mat_SeqAIJ          *a_loc;
1386 
1387   PetscFunctionBegin;
1388   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1389   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1390   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1391 
1392   ptap  = c->ptap;
1393   merge = ptap->merge;
1394 
1395   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1396   /*------------------------------------------*/
1397   /* get data from symbolic products */
1398   coi    = merge->coi; coj = merge->coj;
1399   ierr   = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
1400   ierr   = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
1401   bi     = merge->bi; bj = merge->bj;
1402   owners = merge->rowmap->range;
1403   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1404   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1405 
1406   /* get A_loc by taking all local rows of A */
1407   A_loc = ptap->A_loc;
1408   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1409   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1410   ai    = a_loc->i;
1411   aj    = a_loc->j;
1412 
1413   for (i=0; i<am; i++) {
1414     anz = ai[i+1] - ai[i];
1415     adj = aj + ai[i];
1416     ada = a_loc->a + ai[i];
1417 
1418     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1419     /*-------------------------------------------------------------*/
1420     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1421     pnz = po->i[i+1] - po->i[i];
1422     poJ = po->j + po->i[i];
1423     pA  = po->a + po->i[i];
1424     for (j=0; j<pnz; j++) {
1425       row = poJ[j];
1426       cj  = coj + coi[row];
1427       ca  = coa + coi[row];
1428       /* perform sparse axpy */
1429       nexta  = 0;
1430       valtmp = pA[j];
1431       for (k=0; nexta<anz; k++) {
1432         if (cj[k] == adj[nexta]) {
1433           ca[k] += valtmp*ada[nexta];
1434           nexta++;
1435         }
1436       }
1437       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1438     }
1439 
1440     /* put the value into Cd (diagonal part) */
1441     pnz = pd->i[i+1] - pd->i[i];
1442     pdJ = pd->j + pd->i[i];
1443     pA  = pd->a + pd->i[i];
1444     for (j=0; j<pnz; j++) {
1445       row = pdJ[j];
1446       cj  = bj + bi[row];
1447       ca  = ba + bi[row];
1448       /* perform sparse axpy */
1449       nexta  = 0;
1450       valtmp = pA[j];
1451       for (k=0; nexta<anz; k++) {
1452         if (cj[k] == adj[nexta]) {
1453           ca[k] += valtmp*ada[nexta];
1454           nexta++;
1455         }
1456       }
1457       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1458     }
1459   }
1460 
1461   /* 3) send and recv matrix values coa */
1462   /*------------------------------------*/
1463   buf_ri = merge->buf_ri;
1464   buf_rj = merge->buf_rj;
1465   len_s  = merge->len_s;
1466   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1467   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1468 
1469   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
1470   for (proc=0,k=0; proc<size; proc++) {
1471     if (!len_s[proc]) continue;
1472     i    = merge->owners_co[proc];
1473     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1474     k++;
1475   }
1476   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1477   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1478 
1479   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1480   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1481   ierr = PetscFree(coa);CHKERRQ(ierr);
1482 
1483   /* 4) insert local Cseq and received values into Cmpi */
1484   /*----------------------------------------------------*/
1485   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1486   for (k=0; k<merge->nrecv; k++) {
1487     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1488     nrows       = *(buf_ri_k[k]);
1489     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1490     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1491   }
1492 
1493   for (i=0; i<cm; i++) {
1494     row  = owners[rank] + i; /* global row index of C_seq */
1495     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1496     ba_i = ba + bi[i];
1497     bnz  = bi[i+1] - bi[i];
1498     /* add received vals into ba */
1499     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1500       /* i-th row */
1501       if (i == *nextrow[k]) {
1502         cnz    = *(nextci[k]+1) - *nextci[k];
1503         cj     = buf_rj[k] + *(nextci[k]);
1504         ca     = abuf_r[k] + *(nextci[k]);
1505         nextcj = 0;
1506         for (j=0; nextcj<cnz; j++) {
1507           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1508             ba_i[j] += ca[nextcj++];
1509           }
1510         }
1511         nextrow[k]++; nextci[k]++;
1512         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1513       }
1514     }
1515     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1516   }
1517   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1518   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1519 
1520   ierr = PetscFree(ba);CHKERRQ(ierr);
1521   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1522   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1523   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1524   PetscFunctionReturn(0);
1525 }
1526 
1527 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1528 #undef __FUNCT__
1529 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
1530 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,PetscReal fill,Mat *C)
1531 {
1532   PetscErrorCode      ierr;
1533   Mat                 Cmpi,A_loc,POt,PDt;
1534   Mat_PtAPMPI         *ptap;
1535   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1536   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1537   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1538   PetscInt            nnz;
1539   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1540   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1541   MPI_Comm            comm;
1542   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1543   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1544   PetscInt            len,proc,*dnz,*onz,*owners;
1545   PetscInt            nzi,*bi,*bj;
1546   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1547   MPI_Request         *swaits,*rwaits;
1548   MPI_Status          *sstatus,rstatus;
1549   Mat_Merge_SeqsToMPI *merge;
1550   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1551   PetscReal           afill  =1.0,afill_tmp;
1552   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1553   PetscScalar         *vals;
1554   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1555 
1556   PetscFunctionBegin;
1557   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1558   /* check if matrix local sizes are compatible */
1559   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1560     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);
1561   }
1562 
1563   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1564   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1565 
1566   /* create struct Mat_PtAPMPI and attached it to C later */
1567   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1568 
1569   /* get A_loc by taking all local rows of A */
1570   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1571 
1572   ptap->A_loc = A_loc;
1573   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1574   ai          = a_loc->i;
1575   aj          = a_loc->j;
1576 
1577   /* determine symbolic Co=(p->B)^T*A - send to others */
1578   /*----------------------------------------------------*/
1579   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1580   pdt  = (Mat_SeqAIJ*)PDt->data;
1581   pdti = pdt->i; pdtj = pdt->j;
1582 
1583   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1584   pot  = (Mat_SeqAIJ*)POt->data;
1585   poti = pot->i; potj = pot->j;
1586 
1587   /* then, compute symbolic Co = (p->B)^T*A */
1588   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1589                          >= (num of nonzero rows of C_seq) - pn */
1590   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1591   coi[0] = 0;
1592 
1593   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1594   nnz           = fill*(poti[pon] + ai[am]);
1595   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1596   current_space = free_space;
1597 
1598   /* create and initialize a linked list */
1599   i     = PetscMax(pdt->rmax,pot->rmax);
1600   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1601   if (!Crmax || Crmax > aN) Crmax = aN;
1602   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1603 
1604   for (i=0; i<pon; i++) {
1605     pnz = poti[i+1] - poti[i];
1606     ptJ = potj + poti[i];
1607     for (j=0; j<pnz; j++) {
1608       row  = ptJ[j]; /* row of A_loc == col of Pot */
1609       anz  = ai[row+1] - ai[row];
1610       Jptr = aj + ai[row];
1611       /* add non-zero cols of AP into the sorted linked list lnk */
1612       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1613     }
1614     nnz = lnk[0];
1615 
1616     /* If free space is not available, double the total space in the list */
1617     if (current_space->local_remaining<nnz) {
1618       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1619       nspacedouble++;
1620     }
1621 
1622     /* Copy data into free space, and zero out denserows */
1623     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1624 
1625     current_space->array           += nnz;
1626     current_space->local_used      += nnz;
1627     current_space->local_remaining -= nnz;
1628 
1629     coi[i+1] = coi[i] + nnz;
1630   }
1631 
1632   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1633   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1634 
1635   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1636   if (afill_tmp > afill) afill = afill_tmp;
1637 
1638   /* send j-array (coj) of Co to other processors */
1639   /*----------------------------------------------*/
1640   /* determine row ownership */
1641   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1642   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1643 
1644   merge->rowmap->n  = pn;
1645   merge->rowmap->bs = 1;
1646 
1647   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1648   owners = merge->rowmap->range;
1649 
1650   /* determine the number of messages to send, their lengths */
1651   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1652   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1653   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1654 
1655   len_s        = merge->len_s;
1656   merge->nsend = 0;
1657 
1658   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1659   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1660 
1661   proc = 0;
1662   for (i=0; i<pon; i++) {
1663     while (prmap[i] >= owners[proc+1]) proc++;
1664     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1665     len_s[proc] += coi[i+1] - coi[i];
1666   }
1667 
1668   len          = 0; /* max length of buf_si[] */
1669   owners_co[0] = 0;
1670   for (proc=0; proc<size; proc++) {
1671     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1672     if (len_si[proc]) {
1673       merge->nsend++;
1674       len_si[proc] = 2*(len_si[proc] + 1);
1675       len         += len_si[proc];
1676     }
1677   }
1678 
1679   /* determine the number and length of messages to receive for coi and coj  */
1680   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1681   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1682 
1683   /* post the Irecv and Isend of coj */
1684   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1685   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1686   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1687   for (proc=0, k=0; proc<size; proc++) {
1688     if (!len_s[proc]) continue;
1689     i    = owners_co[proc];
1690     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1691     k++;
1692   }
1693 
1694   /* receives and sends of coj are complete */
1695   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1696   for (i=0; i<merge->nrecv; i++) {
1697     PetscMPIInt icompleted;
1698     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1699   }
1700   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1701   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1702 
1703   /* send and recv coi */
1704   /*-------------------*/
1705   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1706   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1707   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1708   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1709   for (proc=0,k=0; proc<size; proc++) {
1710     if (!len_s[proc]) continue;
1711     /* form outgoing message for i-structure:
1712          buf_si[0]:                 nrows to be sent
1713                [1:nrows]:           row index (global)
1714                [nrows+1:2*nrows+1]: i-structure index
1715     */
1716     /*-------------------------------------------*/
1717     nrows       = len_si[proc]/2 - 1;
1718     buf_si_i    = buf_si + nrows+1;
1719     buf_si[0]   = nrows;
1720     buf_si_i[0] = 0;
1721     nrows       = 0;
1722     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1723       nzi               = coi[i+1] - coi[i];
1724       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1725       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1726       nrows++;
1727     }
1728     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1729     k++;
1730     buf_si += len_si[proc];
1731   }
1732   i = merge->nrecv;
1733   while (i--) {
1734     PetscMPIInt icompleted;
1735     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1736   }
1737   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1738   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1739   ierr = PetscFree(len_si);CHKERRQ(ierr);
1740   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1741   ierr = PetscFree(swaits);CHKERRQ(ierr);
1742   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1743   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1744 
1745   /* compute the local portion of C (mpi mat) */
1746   /*------------------------------------------*/
1747   /* allocate bi array and free space for accumulating nonzero column info */
1748   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1749   bi[0] = 0;
1750 
1751   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1752   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1753   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1754   current_space = free_space;
1755 
1756   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1757   for (k=0; k<merge->nrecv; k++) {
1758     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1759     nrows       = *buf_ri_k[k];
1760     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1761     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1762   }
1763 
1764   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1765   rmax = 0;
1766   for (i=0; i<pn; i++) {
1767     /* add pdt[i,:]*AP into lnk */
1768     pnz = pdti[i+1] - pdti[i];
1769     ptJ = pdtj + pdti[i];
1770     for (j=0; j<pnz; j++) {
1771       row  = ptJ[j];  /* row of AP == col of Pt */
1772       anz  = ai[row+1] - ai[row];
1773       Jptr = aj + ai[row];
1774       /* add non-zero cols of AP into the sorted linked list lnk */
1775       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1776     }
1777 
1778     /* add received col data into lnk */
1779     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1780       if (i == *nextrow[k]) { /* i-th row */
1781         nzi  = *(nextci[k]+1) - *nextci[k];
1782         Jptr = buf_rj[k] + *nextci[k];
1783         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1784         nextrow[k]++; nextci[k]++;
1785       }
1786     }
1787     nnz = lnk[0];
1788 
1789     /* if free space is not available, make more free space */
1790     if (current_space->local_remaining<nnz) {
1791       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1792       nspacedouble++;
1793     }
1794     /* copy data into free space, then initialize lnk */
1795     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1796     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1797 
1798     current_space->array           += nnz;
1799     current_space->local_used      += nnz;
1800     current_space->local_remaining -= nnz;
1801 
1802     bi[i+1] = bi[i] + nnz;
1803     if (nnz > rmax) rmax = nnz;
1804   }
1805   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1806 
1807   ierr      = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1808   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1809   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1810   if (afill_tmp > afill) afill = afill_tmp;
1811   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1812   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1813   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1814 
1815   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1816   /*----------------------------------------------------------------------------------*/
1817   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1818   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1819 
1820   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1821   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1822   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1823   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1824   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1825   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1826   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1827   for (i=0; i<pn; i++) {
1828     row  = i + rstart;
1829     nnz  = bi[i+1] - bi[i];
1830     Jptr = bj + bi[i];
1831     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1832   }
1833   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1834   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1835   ierr = PetscFree(vals);CHKERRQ(ierr);
1836 
1837   merge->bi        = bi;
1838   merge->bj        = bj;
1839   merge->coi       = coi;
1840   merge->coj       = coj;
1841   merge->buf_ri    = buf_ri;
1842   merge->buf_rj    = buf_rj;
1843   merge->owners_co = owners_co;
1844   merge->destroy   = Cmpi->ops->destroy;
1845   merge->duplicate = Cmpi->ops->duplicate;
1846 
1847   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
1848   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1849 
1850   /* attach the supporting struct to Cmpi for reuse */
1851   c = (Mat_MPIAIJ*)Cmpi->data;
1852 
1853   c->ptap     = ptap;
1854   ptap->api   = NULL;
1855   ptap->apj   = NULL;
1856   ptap->merge = merge;
1857   ptap->rmax  = rmax;
1858   ptap->apa   = NULL;
1859 
1860   *C = Cmpi;
1861 #if defined(PETSC_USE_INFO)
1862   if (bi[pn] != 0) {
1863     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1864     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1865   } else {
1866     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1867   }
1868 #endif
1869   PetscFunctionReturn(0);
1870 }
1871