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