xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 0fc8cf3458a76e3df512d381d4ded7ce5e272233)
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_Scalable(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_Scalable"
655 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(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   PetscBool      scalable=PETSC_TRUE,viamatmatmult=PETSC_FALSE;
812 
813   PetscFunctionBegin;
814   ierr = PetscOptionsBool("-mattransposematmult_viamatmatmult","Use R=Pt and C=R*A","",viamatmatmult,&viamatmatmult,NULL);CHKERRQ(ierr);
815   if (viamatmatmult) {
816     Mat         Pt;
817     Mat_PtAPMPI *ptap;
818     Mat_MPIAIJ  *c;
819     if (scall == MAT_INITIAL_MATRIX) {
820       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
821       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
822 
823       c        = (Mat_MPIAIJ*)(*C)->data;
824       ptap     = c->ptap;
825       ptap->Pt = Pt;
826     } else if (scall == MAT_REUSE_MATRIX) {
827       c    = (Mat_MPIAIJ*)(*C)->data;
828       ptap = c->ptap;
829       Pt   = ptap->Pt;
830       ierr = MatTranspose(P,scall,&Pt);CHKERRQ(ierr);
831       ierr = MatMatMult(Pt,A,scall,fill,C);CHKERRQ(ierr);
832     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not supported");
833     PetscFunctionReturn(0);
834   }
835 
836   if (scall == MAT_INITIAL_MATRIX) {
837     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
838     ierr = PetscOptionsBool("-mattransposematmult_scalable","Use a scalable but slower C=Pt*A","",scalable,&scalable,NULL);CHKERRQ(ierr);
839     if  (scalable) {
840       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(P,A,fill,C);CHKERRQ(ierr);
841     } else {
842       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
843     }
844     ierr = PetscOptionsEnd();CHKERRQ(ierr);
845   }
846   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
847   PetscFunctionReturn(0);
848 }
849 
850 #undef __FUNCT__
851 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
852 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
853 {
854   PetscErrorCode      ierr;
855   Mat_Merge_SeqsToMPI *merge;
856   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
857   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
858   Mat_PtAPMPI         *ptap;
859   PetscInt            *adj,*aJ;
860   PetscInt            i,j,k,anz,pnz,row,*cj;
861   MatScalar           *ada,*aval,*ca,valtmp;
862   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
863   MPI_Comm            comm;
864   PetscMPIInt         size,rank,taga,*len_s;
865   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
866   PetscInt            **buf_ri,**buf_rj;
867   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
868   MPI_Request         *s_waits,*r_waits;
869   MPI_Status          *status;
870   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
871   PetscInt            *ai,*aj,*coi,*coj;
872   PetscInt            *poJ,*pdJ;
873   Mat                 A_loc;
874   Mat_SeqAIJ          *a_loc;
875 
876   PetscFunctionBegin;
877   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
878   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
879   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
880 
881   ptap  = c->ptap;
882   merge = ptap->merge;
883 
884   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
885   /*--------------------------------------------------------------*/
886   /* get data from symbolic products */
887   coi  = merge->coi; coj = merge->coj;
888   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
889   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
890 
891   bi     = merge->bi; bj = merge->bj;
892   owners = merge->rowmap->range;
893   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
894   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
895 
896   /* get A_loc by taking all local rows of A */
897   A_loc = ptap->A_loc;
898   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
899   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
900   ai    = a_loc->i;
901   aj    = a_loc->j;
902 
903   ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&aval);CHKERRQ(ierr); /* non-scalable!!! */
904   ierr = PetscMemzero(aval,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
905 
906   for (i=0; i<am; i++) {
907     /* 2-a) put A[i,:] to dense array aval */
908     anz = ai[i+1] - ai[i];
909     adj = aj + ai[i];
910     ada = a_loc->a + ai[i];
911     for (j=0; j<anz; j++) {
912       aval[adj[j]] = ada[j];
913     }
914 
915     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
916     /*--------------------------------------------------------------*/
917     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
918     pnz = po->i[i+1] - po->i[i];
919     poJ = po->j + po->i[i];
920     pA  = po->a + po->i[i];
921     for (j=0; j<pnz; j++) {
922       row = poJ[j];
923       cnz = coi[row+1] - coi[row];
924       cj  = coj + coi[row];
925       ca  = coa + coi[row];
926       /* perform dense axpy */
927       valtmp = pA[j];
928       for (k=0; k<cnz; k++) {
929         ca[k] += valtmp*aval[cj[k]];
930       }
931       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
932     }
933 
934     /* put the value into Cd (diagonal part) */
935     pnz = pd->i[i+1] - pd->i[i];
936     pdJ = pd->j + pd->i[i];
937     pA  = pd->a + pd->i[i];
938     for (j=0; j<pnz; j++) {
939       row = pdJ[j];
940       cnz = bi[row+1] - bi[row];
941       cj  = bj + bi[row];
942       ca  = ba + bi[row];
943       /* perform dense axpy */
944       valtmp = pA[j];
945       for (k=0; k<cnz; k++) {
946         ca[k] += valtmp*aval[cj[k]];
947       }
948       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
949     }
950 
951     /* zero the current row of Pt*A */
952     aJ = aj + ai[i];
953     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
954   }
955 
956   /* 3) send and recv matrix values coa */
957   /*------------------------------------*/
958   buf_ri = merge->buf_ri;
959   buf_rj = merge->buf_rj;
960   len_s  = merge->len_s;
961   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
962   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
963 
964   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
965   for (proc=0,k=0; proc<size; proc++) {
966     if (!len_s[proc]) continue;
967     i    = merge->owners_co[proc];
968     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
969     k++;
970   }
971   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
972   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
973 
974   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
975   ierr = PetscFree(r_waits);CHKERRQ(ierr);
976   ierr = PetscFree(coa);CHKERRQ(ierr);
977 
978   /* 4) insert local Cseq and received values into Cmpi */
979   /*----------------------------------------------------*/
980   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
981   for (k=0; k<merge->nrecv; k++) {
982     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
983     nrows       = *(buf_ri_k[k]);
984     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
985     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
986   }
987 
988   for (i=0; i<cm; i++) {
989     row  = owners[rank] + i; /* global row index of C_seq */
990     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
991     ba_i = ba + bi[i];
992     bnz  = bi[i+1] - bi[i];
993     /* add received vals into ba */
994     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
995       /* i-th row */
996       if (i == *nextrow[k]) {
997         cnz    = *(nextci[k]+1) - *nextci[k];
998         cj     = buf_rj[k] + *(nextci[k]);
999         ca     = abuf_r[k] + *(nextci[k]);
1000         nextcj = 0;
1001         for (j=0; nextcj<cnz; j++) {
1002           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1003             ba_i[j] += ca[nextcj++];
1004           }
1005         }
1006         nextrow[k]++; nextci[k]++;
1007         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1008       }
1009     }
1010     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1011   }
1012   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1013   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1014 
1015   ierr = PetscFree(ba);CHKERRQ(ierr);
1016   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1017   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1018   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1019   ierr = PetscFree(aval);CHKERRQ(ierr);
1020   PetscFunctionReturn(0);
1021 }
1022 
1023 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1024 #undef __FUNCT__
1025 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
1026 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1027 {
1028   PetscErrorCode      ierr;
1029   Mat                 Cmpi,A_loc,POt,PDt;
1030   Mat_PtAPMPI         *ptap;
1031   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1032   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1033   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1034   PetscInt            nnz;
1035   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1036   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1037   PetscBT             lnkbt;
1038   MPI_Comm            comm;
1039   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1040   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1041   PetscInt            len,proc,*dnz,*onz,*owners;
1042   PetscInt            nzi,*bi,*bj;
1043   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1044   MPI_Request         *swaits,*rwaits;
1045   MPI_Status          *sstatus,rstatus;
1046   Mat_Merge_SeqsToMPI *merge;
1047   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1048   PetscReal           afill  =1.0,afill_tmp;
1049   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1050   PetscScalar         *vals;
1051   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1052 
1053   PetscFunctionBegin;
1054   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1055   /* check if matrix local sizes are compatible */
1056   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1057     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);
1058   }
1059 
1060   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1061   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1062 
1063   /* create struct Mat_PtAPMPI and attached it to C later */
1064   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1065 
1066   /* get A_loc by taking all local rows of A */
1067   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1068 
1069   ptap->A_loc = A_loc;
1070 
1071   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1072   ai    = a_loc->i;
1073   aj    = a_loc->j;
1074 
1075   /* determine symbolic Co=(p->B)^T*A - send to others */
1076   /*----------------------------------------------------*/
1077   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1078   pdt  = (Mat_SeqAIJ*)PDt->data;
1079   pdti = pdt->i; pdtj = pdt->j;
1080 
1081   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1082   pot  = (Mat_SeqAIJ*)POt->data;
1083   poti = pot->i; potj = pot->j;
1084 
1085   /* then, compute symbolic Co = (p->B)^T*A */
1086   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1087   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1088   coi[0] = 0;
1089 
1090   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1091   nnz           = fill*(poti[pon] + ai[am]);
1092   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1093   current_space = free_space;
1094 
1095   /* create and initialize a linked list */
1096   i     = PetscMax(pdt->rmax,pot->rmax);
1097   Crmax = i*a_loc->rmax*size;
1098   if (!Crmax || Crmax > aN) Crmax = aN;
1099   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1100 
1101   for (i=0; i<pon; i++) {
1102     pnz = poti[i+1] - poti[i];
1103     ptJ = potj + poti[i];
1104     for (j=0; j<pnz; j++) {
1105       row  = ptJ[j]; /* row of A_loc == col of Pot */
1106       anz  = ai[row+1] - ai[row];
1107       Jptr = aj + ai[row];
1108       /* add non-zero cols of AP into the sorted linked list lnk */
1109       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1110     }
1111     nnz = lnk[0];
1112 
1113     /* If free space is not available, double the total space in the list */
1114     if (current_space->local_remaining<nnz) {
1115       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1116       nspacedouble++;
1117     }
1118 
1119     /* Copy data into free space, and zero out denserows */
1120     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1121 
1122     current_space->array           += nnz;
1123     current_space->local_used      += nnz;
1124     current_space->local_remaining -= nnz;
1125 
1126     coi[i+1] = coi[i] + nnz;
1127   }
1128 
1129   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1130   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1131 
1132   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1133   if (afill_tmp > afill) afill = afill_tmp;
1134 
1135   /* send j-array (coj) of Co to other processors */
1136   /*----------------------------------------------*/
1137   /* determine row ownership */
1138   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1139   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1140 
1141   merge->rowmap->n  = pn;
1142   merge->rowmap->bs = 1;
1143 
1144   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1145   owners = merge->rowmap->range;
1146 
1147   /* determine the number of messages to send, their lengths */
1148   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1149   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1150   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1151 
1152   len_s        = merge->len_s;
1153   merge->nsend = 0;
1154 
1155   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1156   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1157 
1158   proc = 0;
1159   for (i=0; i<pon; i++) {
1160     while (prmap[i] >= owners[proc+1]) proc++;
1161     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1162     len_s[proc] += coi[i+1] - coi[i];
1163   }
1164 
1165   len          = 0; /* max length of buf_si[] */
1166   owners_co[0] = 0;
1167   for (proc=0; proc<size; proc++) {
1168     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1169     if (len_si[proc]) {
1170       merge->nsend++;
1171       len_si[proc] = 2*(len_si[proc] + 1);
1172       len         += len_si[proc];
1173     }
1174   }
1175 
1176   /* determine the number and length of messages to receive for coi and coj  */
1177   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1178   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1179 
1180   /* post the Irecv and Isend of coj */
1181   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1182   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1183   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1184   for (proc=0, k=0; proc<size; proc++) {
1185     if (!len_s[proc]) continue;
1186     i    = owners_co[proc];
1187     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1188     k++;
1189   }
1190 
1191   /* receives and sends of coj are complete */
1192   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1193   for (i=0; i<merge->nrecv; i++) {
1194     PetscMPIInt icompleted;
1195     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1196   }
1197   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1198   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1199 
1200   /* send and recv coi */
1201   /*-------------------*/
1202   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1203   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1204   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1205   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1206   for (proc=0,k=0; proc<size; proc++) {
1207     if (!len_s[proc]) continue;
1208     /* form outgoing message for i-structure:
1209          buf_si[0]:                 nrows to be sent
1210                [1:nrows]:           row index (global)
1211                [nrows+1:2*nrows+1]: i-structure index
1212     */
1213     /*-------------------------------------------*/
1214     nrows       = len_si[proc]/2 - 1;
1215     buf_si_i    = buf_si + nrows+1;
1216     buf_si[0]   = nrows;
1217     buf_si_i[0] = 0;
1218     nrows       = 0;
1219     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1220       nzi               = coi[i+1] - coi[i];
1221       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1222       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1223       nrows++;
1224     }
1225     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1226     k++;
1227     buf_si += len_si[proc];
1228   }
1229   i = merge->nrecv;
1230   while (i--) {
1231     PetscMPIInt icompleted;
1232     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1233   }
1234   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1235   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1236   ierr = PetscFree(len_si);CHKERRQ(ierr);
1237   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1238   ierr = PetscFree(swaits);CHKERRQ(ierr);
1239   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1240   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1241 
1242   /* compute the local portion of C (mpi mat) */
1243   /*------------------------------------------*/
1244   /* allocate bi array and free space for accumulating nonzero column info */
1245   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1246   bi[0] = 0;
1247 
1248   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1249   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1250   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1251   current_space = free_space;
1252 
1253   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1254   for (k=0; k<merge->nrecv; k++) {
1255     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1256     nrows       = *buf_ri_k[k];
1257     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1258     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1259   }
1260 
1261   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1262   rmax = 0;
1263   for (i=0; i<pn; i++) {
1264     /* add pdt[i,:]*AP into lnk */
1265     pnz = pdti[i+1] - pdti[i];
1266     ptJ = pdtj + pdti[i];
1267     for (j=0; j<pnz; j++) {
1268       row  = ptJ[j];  /* row of AP == col of Pt */
1269       anz  = ai[row+1] - ai[row];
1270       Jptr = aj + ai[row];
1271       /* add non-zero cols of AP into the sorted linked list lnk */
1272       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1273     }
1274 
1275     /* add received col data into lnk */
1276     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1277       if (i == *nextrow[k]) { /* i-th row */
1278         nzi  = *(nextci[k]+1) - *nextci[k];
1279         Jptr = buf_rj[k] + *nextci[k];
1280         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1281         nextrow[k]++; nextci[k]++;
1282       }
1283     }
1284     nnz = lnk[0];
1285 
1286     /* if free space is not available, make more free space */
1287     if (current_space->local_remaining<nnz) {
1288       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1289       nspacedouble++;
1290     }
1291     /* copy data into free space, then initialize lnk */
1292     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1293     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1294 
1295     current_space->array           += nnz;
1296     current_space->local_used      += nnz;
1297     current_space->local_remaining -= nnz;
1298 
1299     bi[i+1] = bi[i] + nnz;
1300     if (nnz > rmax) rmax = nnz;
1301   }
1302   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1303 
1304   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1305   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1306 
1307   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1308   if (afill_tmp > afill) afill = afill_tmp;
1309   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
1310   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1311   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1312 
1313   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1314   /*----------------------------------------------------------------------------------*/
1315   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1316   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1317 
1318   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1319   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1320   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1321   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1322   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1323   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1324   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1325   for (i=0; i<pn; i++) {
1326     row  = i + rstart;
1327     nnz  = bi[i+1] - bi[i];
1328     Jptr = bj + bi[i];
1329     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1330   }
1331   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1332   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1333   ierr = PetscFree(vals);CHKERRQ(ierr);
1334 
1335   merge->bi        = bi;
1336   merge->bj        = bj;
1337   merge->coi       = coi;
1338   merge->coj       = coj;
1339   merge->buf_ri    = buf_ri;
1340   merge->buf_rj    = buf_rj;
1341   merge->owners_co = owners_co;
1342   merge->destroy   = Cmpi->ops->destroy;
1343   merge->duplicate = Cmpi->ops->duplicate;
1344 
1345   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1346   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1347 
1348   /* attach the supporting struct to Cmpi for reuse */
1349   c           = (Mat_MPIAIJ*)Cmpi->data;
1350   c->ptap     = ptap;
1351   ptap->api   = NULL;
1352   ptap->apj   = NULL;
1353   ptap->merge = merge;
1354   ptap->rmax  = rmax;
1355 
1356   *C = Cmpi;
1357 #if defined(PETSC_USE_INFO)
1358   if (bi[pn] != 0) {
1359     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1360     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1361   } else {
1362     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1363   }
1364 #endif
1365   PetscFunctionReturn(0);
1366 }
1367 
1368 #undef __FUNCT__
1369 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
1370 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,Mat C)
1371 {
1372   PetscErrorCode      ierr;
1373   Mat_Merge_SeqsToMPI *merge;
1374   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1375   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1376   Mat_PtAPMPI         *ptap;
1377   PetscInt            *adj;
1378   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1379   MatScalar           *ada,*ca,valtmp;
1380   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1381   MPI_Comm            comm;
1382   PetscMPIInt         size,rank,taga,*len_s;
1383   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1384   PetscInt            **buf_ri,**buf_rj;
1385   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1386   MPI_Request         *s_waits,*r_waits;
1387   MPI_Status          *status;
1388   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1389   PetscInt            *ai,*aj,*coi,*coj;
1390   PetscInt            *poJ,*pdJ;
1391   Mat                 A_loc;
1392   Mat_SeqAIJ          *a_loc;
1393 
1394   PetscFunctionBegin;
1395   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1396   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1397   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1398 
1399   ptap  = c->ptap;
1400   merge = ptap->merge;
1401 
1402   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1403   /*------------------------------------------*/
1404   /* get data from symbolic products */
1405   coi    = merge->coi; coj = merge->coj;
1406   ierr   = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
1407   ierr   = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
1408   bi     = merge->bi; bj = merge->bj;
1409   owners = merge->rowmap->range;
1410   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1411   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1412 
1413   /* get A_loc by taking all local rows of A */
1414   A_loc = ptap->A_loc;
1415   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1416   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1417   ai    = a_loc->i;
1418   aj    = a_loc->j;
1419 
1420   for (i=0; i<am; i++) {
1421     anz = ai[i+1] - ai[i];
1422     adj = aj + ai[i];
1423     ada = a_loc->a + ai[i];
1424 
1425     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1426     /*-------------------------------------------------------------*/
1427     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1428     pnz = po->i[i+1] - po->i[i];
1429     poJ = po->j + po->i[i];
1430     pA  = po->a + po->i[i];
1431     for (j=0; j<pnz; j++) {
1432       row = poJ[j];
1433       cj  = coj + coi[row];
1434       ca  = coa + coi[row];
1435       /* perform sparse axpy */
1436       nexta  = 0;
1437       valtmp = pA[j];
1438       for (k=0; nexta<anz; k++) {
1439         if (cj[k] == adj[nexta]) {
1440           ca[k] += valtmp*ada[nexta];
1441           nexta++;
1442         }
1443       }
1444       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1445     }
1446 
1447     /* put the value into Cd (diagonal part) */
1448     pnz = pd->i[i+1] - pd->i[i];
1449     pdJ = pd->j + pd->i[i];
1450     pA  = pd->a + pd->i[i];
1451     for (j=0; j<pnz; j++) {
1452       row = pdJ[j];
1453       cj  = bj + bi[row];
1454       ca  = ba + bi[row];
1455       /* perform sparse axpy */
1456       nexta  = 0;
1457       valtmp = pA[j];
1458       for (k=0; nexta<anz; k++) {
1459         if (cj[k] == adj[nexta]) {
1460           ca[k] += valtmp*ada[nexta];
1461           nexta++;
1462         }
1463       }
1464       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1465     }
1466   }
1467 
1468   /* 3) send and recv matrix values coa */
1469   /*------------------------------------*/
1470   buf_ri = merge->buf_ri;
1471   buf_rj = merge->buf_rj;
1472   len_s  = merge->len_s;
1473   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1474   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1475 
1476   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
1477   for (proc=0,k=0; proc<size; proc++) {
1478     if (!len_s[proc]) continue;
1479     i    = merge->owners_co[proc];
1480     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1481     k++;
1482   }
1483   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1484   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1485 
1486   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1487   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1488   ierr = PetscFree(coa);CHKERRQ(ierr);
1489 
1490   /* 4) insert local Cseq and received values into Cmpi */
1491   /*----------------------------------------------------*/
1492   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1493   for (k=0; k<merge->nrecv; k++) {
1494     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1495     nrows       = *(buf_ri_k[k]);
1496     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1497     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1498   }
1499 
1500   for (i=0; i<cm; i++) {
1501     row  = owners[rank] + i; /* global row index of C_seq */
1502     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1503     ba_i = ba + bi[i];
1504     bnz  = bi[i+1] - bi[i];
1505     /* add received vals into ba */
1506     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1507       /* i-th row */
1508       if (i == *nextrow[k]) {
1509         cnz    = *(nextci[k]+1) - *nextci[k];
1510         cj     = buf_rj[k] + *(nextci[k]);
1511         ca     = abuf_r[k] + *(nextci[k]);
1512         nextcj = 0;
1513         for (j=0; nextcj<cnz; j++) {
1514           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1515             ba_i[j] += ca[nextcj++];
1516           }
1517         }
1518         nextrow[k]++; nextci[k]++;
1519         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1520       }
1521     }
1522     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1523   }
1524   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1525   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1526 
1527   ierr = PetscFree(ba);CHKERRQ(ierr);
1528   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1529   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1530   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1531   PetscFunctionReturn(0);
1532 }
1533 
1534 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1535 #undef __FUNCT__
1536 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
1537 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,PetscReal fill,Mat *C)
1538 {
1539   PetscErrorCode      ierr;
1540   Mat                 Cmpi,A_loc,POt,PDt;
1541   Mat_PtAPMPI         *ptap;
1542   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1543   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1544   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1545   PetscInt            nnz;
1546   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1547   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1548   MPI_Comm            comm;
1549   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1550   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1551   PetscInt            len,proc,*dnz,*onz,*owners;
1552   PetscInt            nzi,*bi,*bj;
1553   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1554   MPI_Request         *swaits,*rwaits;
1555   MPI_Status          *sstatus,rstatus;
1556   Mat_Merge_SeqsToMPI *merge;
1557   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1558   PetscReal           afill  =1.0,afill_tmp;
1559   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1560   PetscScalar         *vals;
1561   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1562 
1563   PetscFunctionBegin;
1564   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1565   /* check if matrix local sizes are compatible */
1566   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1567     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);
1568   }
1569 
1570   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1571   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1572 
1573   /* create struct Mat_PtAPMPI and attached it to C later */
1574   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1575 
1576   /* get A_loc by taking all local rows of A */
1577   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1578 
1579   ptap->A_loc = A_loc;
1580   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1581   ai          = a_loc->i;
1582   aj          = a_loc->j;
1583 
1584   /* determine symbolic Co=(p->B)^T*A - send to others */
1585   /*----------------------------------------------------*/
1586   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1587   pdt  = (Mat_SeqAIJ*)PDt->data;
1588   pdti = pdt->i; pdtj = pdt->j;
1589 
1590   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1591   pot  = (Mat_SeqAIJ*)POt->data;
1592   poti = pot->i; potj = pot->j;
1593 
1594   /* then, compute symbolic Co = (p->B)^T*A */
1595   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1596                          >= (num of nonzero rows of C_seq) - pn */
1597   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1598   coi[0] = 0;
1599 
1600   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1601   nnz           = fill*(poti[pon] + ai[am]);
1602   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1603   current_space = free_space;
1604 
1605   /* create and initialize a linked list */
1606   i     = PetscMax(pdt->rmax,pot->rmax);
1607   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1608   if (!Crmax || Crmax > aN) Crmax = aN;
1609   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1610 
1611   for (i=0; i<pon; i++) {
1612     pnz = poti[i+1] - poti[i];
1613     ptJ = potj + poti[i];
1614     for (j=0; j<pnz; j++) {
1615       row  = ptJ[j]; /* row of A_loc == col of Pot */
1616       anz  = ai[row+1] - ai[row];
1617       Jptr = aj + ai[row];
1618       /* add non-zero cols of AP into the sorted linked list lnk */
1619       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1620     }
1621     nnz = lnk[0];
1622 
1623     /* If free space is not available, double the total space in the list */
1624     if (current_space->local_remaining<nnz) {
1625       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1626       nspacedouble++;
1627     }
1628 
1629     /* Copy data into free space, and zero out denserows */
1630     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1631 
1632     current_space->array           += nnz;
1633     current_space->local_used      += nnz;
1634     current_space->local_remaining -= nnz;
1635 
1636     coi[i+1] = coi[i] + nnz;
1637   }
1638 
1639   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1640   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1641 
1642   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1643   if (afill_tmp > afill) afill = afill_tmp;
1644 
1645   /* send j-array (coj) of Co to other processors */
1646   /*----------------------------------------------*/
1647   /* determine row ownership */
1648   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1649   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1650 
1651   merge->rowmap->n  = pn;
1652   merge->rowmap->bs = 1;
1653 
1654   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1655   owners = merge->rowmap->range;
1656 
1657   /* determine the number of messages to send, their lengths */
1658   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1659   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1660   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1661 
1662   len_s        = merge->len_s;
1663   merge->nsend = 0;
1664 
1665   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1666   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1667 
1668   proc = 0;
1669   for (i=0; i<pon; i++) {
1670     while (prmap[i] >= owners[proc+1]) proc++;
1671     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1672     len_s[proc] += coi[i+1] - coi[i];
1673   }
1674 
1675   len          = 0; /* max length of buf_si[] */
1676   owners_co[0] = 0;
1677   for (proc=0; proc<size; proc++) {
1678     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1679     if (len_si[proc]) {
1680       merge->nsend++;
1681       len_si[proc] = 2*(len_si[proc] + 1);
1682       len         += len_si[proc];
1683     }
1684   }
1685 
1686   /* determine the number and length of messages to receive for coi and coj  */
1687   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1688   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1689 
1690   /* post the Irecv and Isend of coj */
1691   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1692   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1693   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1694   for (proc=0, k=0; proc<size; proc++) {
1695     if (!len_s[proc]) continue;
1696     i    = owners_co[proc];
1697     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1698     k++;
1699   }
1700 
1701   /* receives and sends of coj are complete */
1702   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1703   for (i=0; i<merge->nrecv; i++) {
1704     PetscMPIInt icompleted;
1705     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1706   }
1707   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1708   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1709 
1710   /* send and recv coi */
1711   /*-------------------*/
1712   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1713   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1714   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1715   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1716   for (proc=0,k=0; proc<size; proc++) {
1717     if (!len_s[proc]) continue;
1718     /* form outgoing message for i-structure:
1719          buf_si[0]:                 nrows to be sent
1720                [1:nrows]:           row index (global)
1721                [nrows+1:2*nrows+1]: i-structure index
1722     */
1723     /*-------------------------------------------*/
1724     nrows       = len_si[proc]/2 - 1;
1725     buf_si_i    = buf_si + nrows+1;
1726     buf_si[0]   = nrows;
1727     buf_si_i[0] = 0;
1728     nrows       = 0;
1729     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1730       nzi               = coi[i+1] - coi[i];
1731       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1732       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1733       nrows++;
1734     }
1735     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1736     k++;
1737     buf_si += len_si[proc];
1738   }
1739   i = merge->nrecv;
1740   while (i--) {
1741     PetscMPIInt icompleted;
1742     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1743   }
1744   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1745   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1746   ierr = PetscFree(len_si);CHKERRQ(ierr);
1747   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1748   ierr = PetscFree(swaits);CHKERRQ(ierr);
1749   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1750   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1751 
1752   /* compute the local portion of C (mpi mat) */
1753   /*------------------------------------------*/
1754   /* allocate bi array and free space for accumulating nonzero column info */
1755   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1756   bi[0] = 0;
1757 
1758   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1759   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1760   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1761   current_space = free_space;
1762 
1763   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1764   for (k=0; k<merge->nrecv; k++) {
1765     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1766     nrows       = *buf_ri_k[k];
1767     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1768     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1769   }
1770 
1771   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1772   rmax = 0;
1773   for (i=0; i<pn; i++) {
1774     /* add pdt[i,:]*AP into lnk */
1775     pnz = pdti[i+1] - pdti[i];
1776     ptJ = pdtj + pdti[i];
1777     for (j=0; j<pnz; j++) {
1778       row  = ptJ[j];  /* row of AP == col of Pt */
1779       anz  = ai[row+1] - ai[row];
1780       Jptr = aj + ai[row];
1781       /* add non-zero cols of AP into the sorted linked list lnk */
1782       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1783     }
1784 
1785     /* add received col data into lnk */
1786     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1787       if (i == *nextrow[k]) { /* i-th row */
1788         nzi  = *(nextci[k]+1) - *nextci[k];
1789         Jptr = buf_rj[k] + *nextci[k];
1790         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1791         nextrow[k]++; nextci[k]++;
1792       }
1793     }
1794     nnz = lnk[0];
1795 
1796     /* if free space is not available, make more free space */
1797     if (current_space->local_remaining<nnz) {
1798       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1799       nspacedouble++;
1800     }
1801     /* copy data into free space, then initialize lnk */
1802     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1803     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1804 
1805     current_space->array           += nnz;
1806     current_space->local_used      += nnz;
1807     current_space->local_remaining -= nnz;
1808 
1809     bi[i+1] = bi[i] + nnz;
1810     if (nnz > rmax) rmax = nnz;
1811   }
1812   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1813 
1814   ierr      = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1815   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1816   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1817   if (afill_tmp > afill) afill = afill_tmp;
1818   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1819   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1820   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1821 
1822   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1823   /*----------------------------------------------------------------------------------*/
1824   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1825   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1826 
1827   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1828   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1829   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1830   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1831   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1832   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1833   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1834   for (i=0; i<pn; i++) {
1835     row  = i + rstart;
1836     nnz  = bi[i+1] - bi[i];
1837     Jptr = bj + bi[i];
1838     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1839   }
1840   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1841   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1842   ierr = PetscFree(vals);CHKERRQ(ierr);
1843 
1844   merge->bi        = bi;
1845   merge->bj        = bj;
1846   merge->coi       = coi;
1847   merge->coj       = coj;
1848   merge->buf_ri    = buf_ri;
1849   merge->buf_rj    = buf_rj;
1850   merge->owners_co = owners_co;
1851   merge->destroy   = Cmpi->ops->destroy;
1852   merge->duplicate = Cmpi->ops->duplicate;
1853 
1854   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
1855   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1856 
1857   /* attach the supporting struct to Cmpi for reuse */
1858   c = (Mat_MPIAIJ*)Cmpi->data;
1859 
1860   c->ptap     = ptap;
1861   ptap->api   = NULL;
1862   ptap->apj   = NULL;
1863   ptap->merge = merge;
1864   ptap->rmax  = rmax;
1865   ptap->apa   = NULL;
1866 
1867   *C = Cmpi;
1868 #if defined(PETSC_USE_INFO)
1869   if (bi[pn] != 0) {
1870     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1871     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1872   } else {
1873     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1874   }
1875 #endif
1876   PetscFunctionReturn(0);
1877 }
1878