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