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