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