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