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