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