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