xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 55e7fe800d976e85ed2b5cd8bfdef564daa37bd9)
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;
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   MatType            mtype;
214 
215   PetscFunctionBegin;
216   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
217   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
218 
219   /* create struct Mat_PtAPMPI and attached it to C later */
220   ierr = PetscNew(&ptap);CHKERRQ(ierr);
221 
222   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
223   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
224 
225   /* get P_loc by taking all local rows of P */
226   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
227 
228   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
229   pi_loc = p_loc->i; pj_loc = p_loc->j;
230   if (size > 1) {
231     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
232     pi_oth = p_oth->i; pj_oth = p_oth->j;
233   } else {
234     p_oth = NULL;
235     pi_oth = NULL; pj_oth = NULL;
236   }
237 
238   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
239   /*-------------------------------------------------------------------*/
240   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
241   ptap->api = api;
242   api[0]    = 0;
243 
244   /* create and initialize a linked list */
245   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
246 
247   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
248   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
249   current_space = free_space;
250 
251   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
252   for (i=0; i<am; i++) {
253     /* diagonal portion of A */
254     nzi = adi[i+1] - adi[i];
255     for (j=0; j<nzi; j++) {
256       row  = *adj++;
257       pnz  = pi_loc[row+1] - pi_loc[row];
258       Jptr = pj_loc + pi_loc[row];
259       /* add non-zero cols of P into the sorted linked list lnk */
260       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
261     }
262     /* off-diagonal portion of A */
263     nzi = aoi[i+1] - aoi[i];
264     for (j=0; j<nzi; j++) {
265       row  = *aoj++;
266       pnz  = pi_oth[row+1] - pi_oth[row];
267       Jptr = pj_oth + pi_oth[row];
268       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
269     }
270 
271     apnz     = lnk[0];
272     api[i+1] = api[i] + apnz;
273 
274     /* if free space is not available, double the total space in the list */
275     if (current_space->local_remaining<apnz) {
276       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
277       nspacedouble++;
278     }
279 
280     /* Copy data into free space, then initialize lnk */
281     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
282     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
283 
284     current_space->array           += apnz;
285     current_space->local_used      += apnz;
286     current_space->local_remaining -= apnz;
287   }
288 
289   /* Allocate space for apj, initialize apj, and */
290   /* destroy list of free space and other temporary array(s) */
291   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
292   apj  = ptap->apj;
293   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
294   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
295 
296   /* malloc apa to store dense row A[i,:]*P */
297   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
298 
299   ptap->apa = apa;
300 
301   /* create and assemble symbolic parallel matrix Cmpi */
302   /*----------------------------------------------------*/
303   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
304   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
305   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
306 
307   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
308   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
309   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
310 
311   ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr);
312   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
313   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
314   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
315 
316   ptap->destroy        = Cmpi->ops->destroy;
317   ptap->duplicate      = Cmpi->ops->duplicate;
318   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
319   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
320   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
321 
322   /* attach the supporting struct to Cmpi for reuse */
323   c       = (Mat_MPIAIJ*)Cmpi->data;
324   c->ptap = ptap;
325 
326   *C = Cmpi;
327 
328   /* set MatInfo */
329   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
330   if (afill < 1.0) afill = 1.0;
331   Cmpi->info.mallocs           = nspacedouble;
332   Cmpi->info.fill_ratio_given  = fill;
333   Cmpi->info.fill_ratio_needed = afill;
334 
335 #if defined(PETSC_USE_INFO)
336   if (api[am]) {
337     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
338     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
339   } else {
340     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
341   }
342 #endif
343   PetscFunctionReturn(0);
344 }
345 
346 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
347 {
348   PetscErrorCode ierr;
349 
350   PetscFunctionBegin;
351   if (scall == MAT_INITIAL_MATRIX) {
352     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
353     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
354     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
355   }
356   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
357   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
358   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
359   PetscFunctionReturn(0);
360 }
361 
362 typedef struct {
363   Mat         workB;
364   PetscScalar *rvalues,*svalues;
365   MPI_Request *rwaits,*swaits;
366 } MPIAIJ_MPIDense;
367 
368 PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
369 {
370   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
371   PetscErrorCode  ierr;
372 
373   PetscFunctionBegin;
374   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
375   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
376   ierr = PetscFree(contents);CHKERRQ(ierr);
377   PetscFunctionReturn(0);
378 }
379 
380 #include <petsc/private/vecscatterimpl.h>
381 /*
382     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
383   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
384 
385   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
386 
387   Developer Notes: This directly accesses information inside the VecScatter associated with the matrix-vector product
388    for this matrix. This is not desirable..
389 
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;
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=0;
711   PetscInt           am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20;
712   PetscReal          afill;
713   PetscScalar        *apa;
714   MatType            mtype;
715 
716   PetscFunctionBegin;
717   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
718   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
719 
720   /* create struct Mat_PtAPMPI and attached it to C later */
721   ierr = PetscNew(&ptap);CHKERRQ(ierr);
722 
723   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
724   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
725 
726   /* get P_loc by taking all local rows of P */
727   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
728 
729   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
730   pi_loc = p_loc->i; pj_loc = p_loc->j;
731   if (size > 1) {
732     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
733     pi_oth = p_oth->i; pj_oth = p_oth->j;
734   } else {
735     p_oth  = NULL;
736     pi_oth = NULL; pj_oth = NULL;
737   }
738 
739   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
740   /*-------------------------------------------------------------------*/
741   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
742   ptap->api = api;
743   api[0]    = 0;
744 
745   ierr = PetscLLCondensedCreate_Scalable(lsize,&lnk);CHKERRQ(ierr);
746 
747   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
748   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
749   current_space = free_space;
750   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
751   for (i=0; i<am; i++) {
752     /* diagonal portion of A */
753     nzi = adi[i+1] - adi[i];
754     for (j=0; j<nzi; j++) {
755       row  = *adj++;
756       pnz  = pi_loc[row+1] - pi_loc[row];
757       Jptr = pj_loc + pi_loc[row];
758       /* Expand list if it is not long enough */
759       if (pnz+apnz_max > lsize) {
760         lsize = pnz+apnz_max;
761         ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr);
762       }
763       /* add non-zero cols of P into the sorted linked list lnk */
764       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
765       apnz     = *lnk; /* The first element in the list is the number of items in the list */
766       api[i+1] = api[i] + apnz;
767       if (apnz > apnz_max) apnz_max = apnz;
768     }
769     /* off-diagonal portion of A */
770     nzi = aoi[i+1] - aoi[i];
771     for (j=0; j<nzi; j++) {
772       row  = *aoj++;
773       pnz  = pi_oth[row+1] - pi_oth[row];
774       Jptr = pj_oth + pi_oth[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     apnz     = *lnk;
787     api[i+1] = api[i] + apnz;
788     if (apnz > apnz_max) apnz_max = apnz;
789 
790     /* if free space is not available, double the total space in the list */
791     if (current_space->local_remaining<apnz) {
792       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
793       nspacedouble++;
794     }
795 
796     /* Copy data into free space, then initialize lnk */
797     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
798     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
799 
800     current_space->array           += apnz;
801     current_space->local_used      += apnz;
802     current_space->local_remaining -= apnz;
803   }
804 
805   /* Allocate space for apj, initialize apj, and */
806   /* destroy list of free space and other temporary array(s) */
807   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
808   apj  = ptap->apj;
809   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
810   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
811 
812   /* create and assemble symbolic parallel matrix Cmpi */
813   /*----------------------------------------------------*/
814   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
815   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
816   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
817   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
818   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
819   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
820 
821   /* malloc apa for assembly Cmpi */
822   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
823   ptap->apa = apa;
824 
825   ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr);
826   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
827   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
828   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
829 
830   ptap->destroy             = Cmpi->ops->destroy;
831   ptap->duplicate           = Cmpi->ops->duplicate;
832   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
833   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
834   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
835 
836   /* attach the supporting struct to Cmpi for reuse */
837   c       = (Mat_MPIAIJ*)Cmpi->data;
838   c->ptap = ptap;
839   *C = Cmpi;
840 
841   /* set MatInfo */
842   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
843   if (afill < 1.0) afill = 1.0;
844   Cmpi->info.mallocs           = nspacedouble;
845   Cmpi->info.fill_ratio_given  = fill;
846   Cmpi->info.fill_ratio_needed = afill;
847 
848 #if defined(PETSC_USE_INFO)
849   if (api[am]) {
850     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
851     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
852   } else {
853     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
854   }
855 #endif
856   PetscFunctionReturn(0);
857 }
858 
859 /* This function is needed for the seqMPI matrix-matrix multiplication.  */
860 /* Three input arrays are merged to one output array. The size of the    */
861 /* output array is also output. Duplicate entries only show up once.     */
862 static void Merge3SortedArrays(PetscInt  size1, PetscInt *in1,
863                                PetscInt  size2, PetscInt *in2,
864                                PetscInt  size3, PetscInt *in3,
865                                PetscInt *size4, PetscInt *out)
866 {
867   int i = 0, j = 0, k = 0, l = 0;
868 
869   /* Traverse all three arrays */
870   while (i<size1 && j<size2 && k<size3) {
871     if (in1[i] < in2[j] && in1[i] < in3[k]) {
872       out[l++] = in1[i++];
873     }
874     else if(in2[j] < in1[i] && in2[j] < in3[k]) {
875       out[l++] = in2[j++];
876     }
877     else if(in3[k] < in1[i] && in3[k] < in2[j]) {
878       out[l++] = in3[k++];
879     }
880     else if(in1[i] == in2[j] && in1[i] < in3[k]) {
881       out[l++] = in1[i];
882       i++, j++;
883     }
884     else if(in1[i] == in3[k] && in1[i] < in2[j]) {
885       out[l++] = in1[i];
886       i++, k++;
887     }
888     else if(in3[k] == in2[j] && in2[j] < in1[i])  {
889       out[l++] = in2[j];
890       k++, j++;
891     }
892     else if(in1[i] == in2[j] && in1[i] == in3[k]) {
893       out[l++] = in1[i];
894       i++, j++, k++;
895     }
896   }
897 
898   /* Traverse two remaining arrays */
899   while (i<size1 && j<size2) {
900     if (in1[i] < in2[j]) {
901       out[l++] = in1[i++];
902     }
903     else if(in1[i] > in2[j]) {
904       out[l++] = in2[j++];
905     }
906     else {
907       out[l++] = in1[i];
908       i++, j++;
909     }
910   }
911 
912   while (i<size1 && k<size3) {
913     if (in1[i] < in3[k]) {
914       out[l++] = in1[i++];
915     }
916     else if(in1[i] > in3[k]) {
917       out[l++] = in3[k++];
918     }
919     else {
920       out[l++] = in1[i];
921       i++, k++;
922     }
923   }
924 
925   while (k<size3 && j<size2)  {
926     if (in3[k] < in2[j]) {
927       out[l++] = in3[k++];
928     }
929     else if(in3[k] > in2[j]) {
930       out[l++] = in2[j++];
931     }
932     else {
933       out[l++] = in3[k];
934       k++, j++;
935     }
936   }
937 
938   /* Traverse one remaining array */
939   while (i<size1) out[l++] = in1[i++];
940   while (j<size2) out[l++] = in2[j++];
941   while (k<size3) out[l++] = in3[k++];
942 
943   *size4 = l;
944 }
945 
946 /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and  */
947 /* adds up the products. Two of these three multiplications are performed with existing (sequential)      */
948 /* matrix-matrix multiplications.  */
949 #undef __FUNCT__
950 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI"
951 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat *C)
952 {
953   PetscErrorCode     ierr;
954   MPI_Comm           comm;
955   PetscMPIInt        size;
956   Mat                Cmpi;
957   Mat_PtAPMPI        *ptap;
958   PetscFreeSpaceList free_space_diag=NULL, current_space=NULL;
959   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data;
960   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc;
961   Mat_MPIAIJ         *p        =(Mat_MPIAIJ*)P->data;
962   Mat_MPIAIJ         *c;
963   Mat_SeqAIJ         *adpd_seq, *p_off, *aopoth_seq;
964   PetscInt           adponz, adpdnz;
965   PetscInt           *pi_loc,*dnz,*onz;
966   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart;
967   PetscInt           *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi,
968                      *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj;
969   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend;
970   PetscBT            lnkbt;
971   PetscScalar        *apa;
972   PetscReal          afill;
973   PetscMPIInt        rank;
974   Mat                adpd, aopoth;
975   MatType            mtype;
976 
977   PetscFunctionBegin;
978   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
979   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
980   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
981   ierr = MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend); CHKERRQ(ierr);
982 
983   /* create struct Mat_PtAPMPI and attached it to C later */
984   ierr = PetscNew(&ptap);CHKERRQ(ierr);
985 
986   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
987   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
988 
989   /* get P_loc by taking all local rows of P */
990   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
991 
992 
993   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
994   pi_loc = p_loc->i;
995 
996   /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */
997   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
998   ierr      = PetscMalloc1(am+2,&adpoi);CHKERRQ(ierr);
999 
1000   adpoi[0]    = 0;
1001   ptap->api = api;
1002   api[0] = 0;
1003 
1004   /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */
1005   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
1006   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
1007 
1008   /* Symbolic calc of A_loc_diag * P_loc_diag */
1009   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->A, p->A, fill, &adpd);CHKERRQ(ierr);
1010   adpd_seq = (Mat_SeqAIJ*)((adpd)->data);
1011   adpdi = adpd_seq->i; adpdj = adpd_seq->j;
1012   p_off = (Mat_SeqAIJ*)((p->B)->data);
1013   poff_i = p_off->i; poff_j = p_off->j;
1014 
1015   /* j_temp stores indices of a result row before they are added to the linked list */
1016   ierr = PetscMalloc1(pN+2,&j_temp);CHKERRQ(ierr);
1017 
1018 
1019   /* Symbolic calc of the A_diag * p_loc_off */
1020   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
1021   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);CHKERRQ(ierr);
1022   current_space = free_space_diag;
1023 
1024   for (i=0; i<am; i++) {
1025     /* A_diag * P_loc_off */
1026     nzi = adi[i+1] - adi[i];
1027     for (j=0; j<nzi; j++) {
1028       row  = *adj++;
1029       pnz  = poff_i[row+1] - poff_i[row];
1030       Jptr = poff_j + poff_i[row];
1031       for(i1 = 0; i1 < pnz; i1++) {
1032         j_temp[i1] = p->garray[Jptr[i1]];
1033       }
1034       /* add non-zero cols of P into the sorted linked list lnk */
1035       ierr = PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);CHKERRQ(ierr);
1036     }
1037 
1038     adponz     = lnk[0];
1039     adpoi[i+1] = adpoi[i] + adponz;
1040 
1041     /* if free space is not available, double the total space in the list */
1042     if (current_space->local_remaining<adponz) {
1043       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1044       nspacedouble++;
1045     }
1046 
1047     /* Copy data into free space, then initialize lnk */
1048     ierr = PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1049 
1050     current_space->array           += adponz;
1051     current_space->local_used      += adponz;
1052     current_space->local_remaining -= adponz;
1053   }
1054 
1055   /* Symbolic calc of A_off * P_oth */
1056   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, &aopoth);CHKERRQ(ierr);
1057   aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data);
1058   aopothi = aopoth_seq->i; aopothj = aopoth_seq->j;
1059 
1060   /* Allocate space for apj, adpj, aopj, ... */
1061   /* destroy lists of free space and other temporary array(s) */
1062 
1063   ierr = PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);CHKERRQ(ierr);
1064   ierr = PetscMalloc1(adpoi[am]+2, &adpoj);CHKERRQ(ierr);
1065 
1066   /* Copy from linked list to j-array */
1067   ierr = PetscFreeSpaceContiguous(&free_space_diag,adpoj);CHKERRQ(ierr);
1068   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1069 
1070   adpoJ = adpoj;
1071   adpdJ = adpdj;
1072   aopJ = aopothj;
1073   apj  = ptap->apj;
1074   apJ = apj; /* still empty */
1075 
1076   /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1077   /* A_diag * P_loc_diag to get A*P */
1078   for (i = 0; i < am; i++) {
1079     aopnz  =  aopothi[i+1] -  aopothi[i];
1080     adponz = adpoi[i+1] - adpoi[i];
1081     adpdnz = adpdi[i+1] - adpdi[i];
1082 
1083     /* Correct indices from A_diag*P_diag */
1084     for(i1 = 0; i1 < adpdnz; i1++) {
1085       adpdJ[i1] += p_colstart;
1086     }
1087     /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1088     Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1089     ierr = MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz); CHKERRQ(ierr);
1090 
1091     aopJ += aopnz;
1092     adpoJ += adponz;
1093     adpdJ += adpdnz;
1094     apJ += apnz;
1095     api[i+1] = api[i] + apnz;
1096   }
1097 
1098   /* malloc apa to store dense row A[i,:]*P */
1099   ierr = PetscCalloc1(pN+2,&apa);CHKERRQ(ierr);
1100 
1101   ptap->apa = apa;
1102   /* create and assemble symbolic parallel matrix Cmpi */
1103   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1104   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1105   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
1106   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1107   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
1108   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1109 
1110 
1111   ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr);
1112   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1113   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1114   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1115 
1116 
1117   ptap->destroy        = Cmpi->ops->destroy;
1118   ptap->duplicate      = Cmpi->ops->duplicate;
1119   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1120   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
1121   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
1122 
1123   /* attach the supporting struct to Cmpi for reuse */
1124   c       = (Mat_MPIAIJ*)Cmpi->data;
1125   c->ptap = ptap;
1126   *C = Cmpi;
1127 
1128   /* set MatInfo */
1129   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
1130   if (afill < 1.0) afill = 1.0;
1131   Cmpi->info.mallocs           = nspacedouble;
1132   Cmpi->info.fill_ratio_given  = fill;
1133   Cmpi->info.fill_ratio_needed = afill;
1134 
1135 #if defined(PETSC_USE_INFO)
1136   if (api[am]) {
1137     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
1138     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
1139   } else {
1140     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1141   }
1142 #endif
1143 
1144   ierr = MatDestroy(&aopoth);CHKERRQ(ierr);
1145   ierr = MatDestroy(&adpd);CHKERRQ(ierr);
1146   ierr = PetscFree(j_temp);CHKERRQ(ierr);
1147   ierr = PetscFree(adpoj);CHKERRQ(ierr);
1148   ierr = PetscFree(adpoi);CHKERRQ(ierr);
1149   PetscFunctionReturn(0);
1150 }
1151 
1152 
1153 /*-------------------------------------------------------------------------*/
1154 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
1155 {
1156   PetscErrorCode ierr;
1157   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
1158   PetscInt       aN=A->cmap->N,alg=1; /* set default algorithm */
1159   PetscBool      flg;
1160 
1161   PetscFunctionBegin;
1162   if (scall == MAT_INITIAL_MATRIX) {
1163     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
1164     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
1165     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[1],&alg,&flg);CHKERRQ(ierr);
1166     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1167 
1168     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
1169     switch (alg) {
1170     case 1:
1171       if (!flg && aN > 100000) { /* may switch to scalable algorithm as default */
1172         MatInfo     Ainfo,Pinfo;
1173         PetscInt    nz_local;
1174         PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
1175         MPI_Comm    comm;
1176 
1177         ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
1178         ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
1179         nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); /* estimated local nonzero entries */
1180 
1181         if (aN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
1182         ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1183         ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
1184 
1185         if (alg_scalable) {
1186           alg  = 0; /* scalable algorithm would slower than nonscalable algorithm */
1187           ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
1188           break;
1189         }
1190       }
1191       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
1192       break;
1193     case 2:
1194     {
1195       Mat         Pt;
1196       Mat_PtAPMPI *ptap;
1197       Mat_MPIAIJ  *c;
1198       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
1199       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
1200       c        = (Mat_MPIAIJ*)(*C)->data;
1201       ptap     = c->ptap;
1202       ptap->Pt = Pt;
1203       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
1204       PetscFunctionReturn(0);
1205     }
1206       break;
1207     default: /* scalable algorithm */
1208       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
1209       break;
1210     }
1211     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
1212   }
1213   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
1214   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
1215   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
1216   PetscFunctionReturn(0);
1217 }
1218 
1219 /* This routine only works when scall=MAT_REUSE_MATRIX! */
1220 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
1221 {
1222   PetscErrorCode ierr;
1223   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
1224   Mat_PtAPMPI    *ptap= c->ptap;
1225   Mat            Pt=ptap->Pt;
1226 
1227   PetscFunctionBegin;
1228   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
1229   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
1230   PetscFunctionReturn(0);
1231 }
1232 
1233 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat,MatDuplicateOption,Mat*);
1234 
1235 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1236 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1237 {
1238   PetscErrorCode      ierr;
1239   Mat_PtAPMPI         *ptap;
1240   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*c;
1241   MPI_Comm            comm;
1242   PetscMPIInt         size,rank;
1243   Mat                 Cmpi;
1244   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1245   PetscInt            pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n;
1246   PetscInt            *lnk,i,k,nsend;
1247   PetscBT             lnkbt;
1248   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
1249   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1250   PetscInt            len,proc,*dnz,*onz,*owners,nzi;
1251   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1252   MPI_Request         *swaits,*rwaits;
1253   MPI_Status          *sstatus,rstatus;
1254   PetscLayout         rowmap;
1255   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
1256   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
1257   PetscInt            *Jptr,*prmap=p->garray,con,j,Crmax;
1258   Mat_SeqAIJ          *a_loc,*c_loc,*c_oth;
1259   PetscTable          ta;
1260   MatType             mtype;
1261 
1262   PetscFunctionBegin;
1263   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1264   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1265   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1266 
1267   /* create symbolic parallel matrix Cmpi */
1268   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1269   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1270   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
1271 
1272   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1273 
1274   /* create struct Mat_PtAPMPI and attached it to C later */
1275   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1276   ptap->reuse = MAT_INITIAL_MATRIX;
1277 
1278   /* (0) compute Rd = Pd^T, Ro = Po^T  */
1279   /* --------------------------------- */
1280   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1281   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1282 
1283   /* (1) compute symbolic A_loc */
1284   /* ---------------------------*/
1285   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);CHKERRQ(ierr);
1286 
1287   /* (2-1) compute symbolic C_oth = Ro*A_loc  */
1288   /* ------------------------------------ */
1289   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
1290 
1291   /* (3) send coj of C_oth to other processors  */
1292   /* ------------------------------------------ */
1293   /* determine row ownership */
1294   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
1295   rowmap->n  = pn;
1296   rowmap->bs = 1;
1297   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
1298   owners = rowmap->range;
1299 
1300   /* determine the number of messages to send, their lengths */
1301   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
1302   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1303   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1304 
1305   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1306   coi   = c_oth->i; coj = c_oth->j;
1307   con   = ptap->C_oth->rmap->n;
1308   proc  = 0;
1309   for (i=0; i<con; i++) {
1310     while (prmap[i] >= owners[proc+1]) proc++;
1311     len_si[proc]++;               /* num of rows in Co(=Pt*A) to be sent to [proc] */
1312     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1313   }
1314 
1315   len          = 0; /* max length of buf_si[], see (4) */
1316   owners_co[0] = 0;
1317   nsend        = 0;
1318   for (proc=0; proc<size; proc++) {
1319     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1320     if (len_s[proc]) {
1321       nsend++;
1322       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1323       len         += len_si[proc];
1324     }
1325   }
1326 
1327   /* determine the number and length of messages to receive for coi and coj  */
1328   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
1329   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
1330 
1331   /* post the Irecv and Isend of coj */
1332   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1333   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1334   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
1335   for (proc=0, k=0; proc<size; proc++) {
1336     if (!len_s[proc]) continue;
1337     i    = owners_co[proc];
1338     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1339     k++;
1340   }
1341 
1342   /* (2-2) compute symbolic C_loc = Rd*A_loc */
1343   /* ---------------------------------------- */
1344   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
1345   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
1346 
1347   /* receives coj are complete */
1348   for (i=0; i<nrecv; i++) {
1349     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1350   }
1351   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1352   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
1353 
1354   /* add received column indices into ta to update Crmax */
1355   a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data;
1356 
1357   /* create and initialize a linked list */
1358   ierr = PetscTableCreate(an,aN,&ta);CHKERRQ(ierr); /* for compute Crmax */
1359   MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta);
1360 
1361   for (k=0; k<nrecv; k++) {/* k-th received message */
1362     Jptr = buf_rj[k];
1363     for (j=0; j<len_r[k]; j++) {
1364       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
1365     }
1366   }
1367   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
1368   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
1369 
1370   /* (4) send and recv coi */
1371   /*-----------------------*/
1372   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1373   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1374   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1375   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1376   for (proc=0,k=0; proc<size; proc++) {
1377     if (!len_s[proc]) continue;
1378     /* form outgoing message for i-structure:
1379          buf_si[0]:                 nrows to be sent
1380                [1:nrows]:           row index (global)
1381                [nrows+1:2*nrows+1]: i-structure index
1382     */
1383     /*-------------------------------------------*/
1384     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1385     buf_si_i    = buf_si + nrows+1;
1386     buf_si[0]   = nrows;
1387     buf_si_i[0] = 0;
1388     nrows       = 0;
1389     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1390       nzi = coi[i+1] - coi[i];
1391       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1392       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1393       nrows++;
1394     }
1395     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1396     k++;
1397     buf_si += len_si[proc];
1398   }
1399   for (i=0; i<nrecv; i++) {
1400     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1401   }
1402   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1403   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
1404 
1405   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
1406   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1407   ierr = PetscFree(swaits);CHKERRQ(ierr);
1408   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1409 
1410   /* (5) compute the local portion of Cmpi      */
1411   /* ------------------------------------------ */
1412   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
1413   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
1414   current_space = free_space;
1415 
1416   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
1417   for (k=0; k<nrecv; k++) {
1418     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1419     nrows       = *buf_ri_k[k];
1420     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1421     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1422   }
1423 
1424   ierr = MatPreallocateInitialize(comm,pn,an,dnz,onz);CHKERRQ(ierr);
1425   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1426   for (i=0; i<pn; i++) {
1427     /* add C_loc into Cmpi */
1428     nzi  = c_loc->i[i+1] - c_loc->i[i];
1429     Jptr = c_loc->j + c_loc->i[i];
1430     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1431 
1432     /* add received col data into lnk */
1433     for (k=0; k<nrecv; k++) { /* k-th received message */
1434       if (i == *nextrow[k]) { /* i-th row */
1435         nzi  = *(nextci[k]+1) - *nextci[k];
1436         Jptr = buf_rj[k] + *nextci[k];
1437         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1438         nextrow[k]++; nextci[k]++;
1439       }
1440     }
1441     nzi = lnk[0];
1442 
1443     /* copy data into free space, then initialize lnk */
1444     ierr = PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1445     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
1446   }
1447   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1448   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1449   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
1450 
1451   /* local sizes and preallocation */
1452   ierr = MatSetSizes(Cmpi,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1453   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
1454   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1455   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1456 
1457   /* members in merge */
1458   ierr = PetscFree(id_r);CHKERRQ(ierr);
1459   ierr = PetscFree(len_r);CHKERRQ(ierr);
1460   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
1461   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
1462   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
1463   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
1464   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
1465 
1466   /* attach the supporting struct to Cmpi for reuse */
1467   c = (Mat_MPIAIJ*)Cmpi->data;
1468   c->ptap         = ptap;
1469   ptap->duplicate = Cmpi->ops->duplicate;
1470   ptap->destroy   = Cmpi->ops->destroy;
1471 
1472   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1473   Cmpi->assembled        = PETSC_FALSE;
1474   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1475   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
1476   *C                     = Cmpi;
1477   PetscFunctionReturn(0);
1478 }
1479 
1480 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
1481 {
1482   PetscErrorCode    ierr;
1483   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1484   Mat_SeqAIJ        *c_seq;
1485   Mat_PtAPMPI       *ptap = c->ptap;
1486   Mat               A_loc,C_loc,C_oth;
1487   PetscInt          i,rstart,rend,cm,ncols,row;
1488   const PetscInt    *cols;
1489   const PetscScalar *vals;
1490 
1491   PetscFunctionBegin;
1492   ierr = MatZeroEntries(C);CHKERRQ(ierr);
1493 
1494   if (ptap->reuse == MAT_REUSE_MATRIX) {
1495     /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1496     /* 1) get R = Pd^T, Ro = Po^T */
1497     /*----------------------------*/
1498     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1499     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1500 
1501     /* 2) compute numeric A_loc */
1502     /*--------------------------*/
1503     ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);CHKERRQ(ierr);
1504   }
1505 
1506   /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
1507   A_loc = ptap->A_loc;
1508   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);CHKERRQ(ierr);
1509   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);CHKERRQ(ierr);
1510   C_loc = ptap->C_loc;
1511   C_oth = ptap->C_oth;
1512 
1513   /* add C_loc and Co to to C */
1514   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
1515 
1516   /* C_loc -> C */
1517   cm    = C_loc->rmap->N;
1518   c_seq = (Mat_SeqAIJ*)C_loc->data;
1519   cols = c_seq->j;
1520   vals = c_seq->a;
1521   for (i=0; i<cm; i++) {
1522     ncols = c_seq->i[i+1] - c_seq->i[i];
1523     row = rstart + i;
1524     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1525     cols += ncols; vals += ncols;
1526   }
1527 
1528   /* Co -> C, off-processor part */
1529   cm    = C_oth->rmap->N;
1530   c_seq = (Mat_SeqAIJ*)C_oth->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 = p->garray[i];
1536     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1537     cols += ncols; vals += ncols;
1538   }
1539   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1540   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1541 
1542   ptap->reuse = MAT_REUSE_MATRIX;
1543   PetscFunctionReturn(0);
1544 }
1545 
1546 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1547 {
1548   PetscErrorCode      ierr;
1549   Mat_Merge_SeqsToMPI *merge;
1550   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1551   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1552   Mat_PtAPMPI         *ptap;
1553   PetscInt            *adj;
1554   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1555   MatScalar           *ada,*ca,valtmp;
1556   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1557   MPI_Comm            comm;
1558   PetscMPIInt         size,rank,taga,*len_s;
1559   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1560   PetscInt            **buf_ri,**buf_rj;
1561   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1562   MPI_Request         *s_waits,*r_waits;
1563   MPI_Status          *status;
1564   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1565   PetscInt            *ai,*aj,*coi,*coj,*poJ,*pdJ;
1566   Mat                 A_loc;
1567   Mat_SeqAIJ          *a_loc;
1568 
1569   PetscFunctionBegin;
1570   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1571   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1572   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1573 
1574   ptap  = c->ptap;
1575   merge = ptap->merge;
1576 
1577   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1578   /*------------------------------------------*/
1579   /* get data from symbolic products */
1580   coi    = merge->coi; coj = merge->coj;
1581   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1582   bi     = merge->bi; bj = merge->bj;
1583   owners = merge->rowmap->range;
1584   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1585 
1586   /* get A_loc by taking all local rows of A */
1587   A_loc = ptap->A_loc;
1588   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1589   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1590   ai    = a_loc->i;
1591   aj    = a_loc->j;
1592 
1593   for (i=0; i<am; i++) {
1594     anz = ai[i+1] - ai[i];
1595     adj = aj + ai[i];
1596     ada = a_loc->a + ai[i];
1597 
1598     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1599     /*-------------------------------------------------------------*/
1600     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1601     pnz = po->i[i+1] - po->i[i];
1602     poJ = po->j + po->i[i];
1603     pA  = po->a + po->i[i];
1604     for (j=0; j<pnz; j++) {
1605       row = poJ[j];
1606       cj  = coj + coi[row];
1607       ca  = coa + coi[row];
1608       /* perform sparse axpy */
1609       nexta  = 0;
1610       valtmp = pA[j];
1611       for (k=0; nexta<anz; k++) {
1612         if (cj[k] == adj[nexta]) {
1613           ca[k] += valtmp*ada[nexta];
1614           nexta++;
1615         }
1616       }
1617       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1618     }
1619 
1620     /* put the value into Cd (diagonal part) */
1621     pnz = pd->i[i+1] - pd->i[i];
1622     pdJ = pd->j + pd->i[i];
1623     pA  = pd->a + pd->i[i];
1624     for (j=0; j<pnz; j++) {
1625       row = pdJ[j];
1626       cj  = bj + bi[row];
1627       ca  = ba + bi[row];
1628       /* perform sparse axpy */
1629       nexta  = 0;
1630       valtmp = pA[j];
1631       for (k=0; nexta<anz; k++) {
1632         if (cj[k] == adj[nexta]) {
1633           ca[k] += valtmp*ada[nexta];
1634           nexta++;
1635         }
1636       }
1637       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1638     }
1639   }
1640 
1641   /* 3) send and recv matrix values coa */
1642   /*------------------------------------*/
1643   buf_ri = merge->buf_ri;
1644   buf_rj = merge->buf_rj;
1645   len_s  = merge->len_s;
1646   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1647   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1648 
1649   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1650   for (proc=0,k=0; proc<size; proc++) {
1651     if (!len_s[proc]) continue;
1652     i    = merge->owners_co[proc];
1653     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1654     k++;
1655   }
1656   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1657   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1658 
1659   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1660   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1661   ierr = PetscFree(coa);CHKERRQ(ierr);
1662 
1663   /* 4) insert local Cseq and received values into Cmpi */
1664   /*----------------------------------------------------*/
1665   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1666   for (k=0; k<merge->nrecv; k++) {
1667     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1668     nrows       = *(buf_ri_k[k]);
1669     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1670     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1671   }
1672 
1673   for (i=0; i<cm; i++) {
1674     row  = owners[rank] + i; /* global row index of C_seq */
1675     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1676     ba_i = ba + bi[i];
1677     bnz  = bi[i+1] - bi[i];
1678     /* add received vals into ba */
1679     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1680       /* i-th row */
1681       if (i == *nextrow[k]) {
1682         cnz    = *(nextci[k]+1) - *nextci[k];
1683         cj     = buf_rj[k] + *(nextci[k]);
1684         ca     = abuf_r[k] + *(nextci[k]);
1685         nextcj = 0;
1686         for (j=0; nextcj<cnz; j++) {
1687           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1688             ba_i[j] += ca[nextcj++];
1689           }
1690         }
1691         nextrow[k]++; nextci[k]++;
1692         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1693       }
1694     }
1695     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1696   }
1697   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1698   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1699 
1700   ierr = PetscFree(ba);CHKERRQ(ierr);
1701   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1702   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1703   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1704   PetscFunctionReturn(0);
1705 }
1706 
1707 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1708 {
1709   PetscErrorCode      ierr;
1710   Mat                 Cmpi,A_loc,POt,PDt;
1711   Mat_PtAPMPI         *ptap;
1712   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1713   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1714   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1715   PetscInt            nnz;
1716   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1717   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1718   MPI_Comm            comm;
1719   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1720   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1721   PetscInt            len,proc,*dnz,*onz,*owners;
1722   PetscInt            nzi,*bi,*bj;
1723   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1724   MPI_Request         *swaits,*rwaits;
1725   MPI_Status          *sstatus,rstatus;
1726   Mat_Merge_SeqsToMPI *merge;
1727   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1728   PetscReal           afill  =1.0,afill_tmp;
1729   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1730   PetscScalar         *vals;
1731   Mat_SeqAIJ          *a_loc,*pdt,*pot;
1732   PetscTable          ta;
1733   MatType             mtype;
1734 
1735   PetscFunctionBegin;
1736   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1737   /* check if matrix local sizes are compatible */
1738   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);
1739 
1740   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1741   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1742 
1743   /* create struct Mat_PtAPMPI and attached it to C later */
1744   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1745 
1746   /* get A_loc by taking all local rows of A */
1747   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1748 
1749   ptap->A_loc = A_loc;
1750   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1751   ai          = a_loc->i;
1752   aj          = a_loc->j;
1753 
1754   /* determine symbolic Co=(p->B)^T*A - send to others */
1755   /*----------------------------------------------------*/
1756   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1757   pdt  = (Mat_SeqAIJ*)PDt->data;
1758   pdti = pdt->i; pdtj = pdt->j;
1759 
1760   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1761   pot  = (Mat_SeqAIJ*)POt->data;
1762   poti = pot->i; potj = pot->j;
1763 
1764   /* then, compute symbolic Co = (p->B)^T*A */
1765   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1766                          >= (num of nonzero rows of C_seq) - pn */
1767   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1768   coi[0] = 0;
1769 
1770   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1771   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1772   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1773   current_space = free_space;
1774 
1775   /* create and initialize a linked list */
1776   ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr);
1777   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
1778   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
1779 
1780   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1781 
1782   for (i=0; i<pon; i++) {
1783     pnz = poti[i+1] - poti[i];
1784     ptJ = potj + poti[i];
1785     for (j=0; j<pnz; j++) {
1786       row  = ptJ[j]; /* row of A_loc == col of Pot */
1787       anz  = ai[row+1] - ai[row];
1788       Jptr = aj + ai[row];
1789       /* add non-zero cols of AP into the sorted linked list lnk */
1790       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1791     }
1792     nnz = lnk[0];
1793 
1794     /* If free space is not available, double the total space in the list */
1795     if (current_space->local_remaining<nnz) {
1796       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1797       nspacedouble++;
1798     }
1799 
1800     /* Copy data into free space, and zero out denserows */
1801     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1802 
1803     current_space->array           += nnz;
1804     current_space->local_used      += nnz;
1805     current_space->local_remaining -= nnz;
1806 
1807     coi[i+1] = coi[i] + nnz;
1808   }
1809 
1810   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1811   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1812   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
1813 
1814   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1815   if (afill_tmp > afill) afill = afill_tmp;
1816 
1817   /* send j-array (coj) of Co to other processors */
1818   /*----------------------------------------------*/
1819   /* determine row ownership */
1820   ierr = PetscNew(&merge);CHKERRQ(ierr);
1821   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1822 
1823   merge->rowmap->n  = pn;
1824   merge->rowmap->bs = 1;
1825 
1826   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1827   owners = merge->rowmap->range;
1828 
1829   /* determine the number of messages to send, their lengths */
1830   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1831   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
1832 
1833   len_s        = merge->len_s;
1834   merge->nsend = 0;
1835 
1836   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1837   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1838 
1839   proc = 0;
1840   for (i=0; i<pon; i++) {
1841     while (prmap[i] >= owners[proc+1]) proc++;
1842     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1843     len_s[proc] += coi[i+1] - coi[i];
1844   }
1845 
1846   len          = 0; /* max length of buf_si[] */
1847   owners_co[0] = 0;
1848   for (proc=0; proc<size; proc++) {
1849     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1850     if (len_si[proc]) {
1851       merge->nsend++;
1852       len_si[proc] = 2*(len_si[proc] + 1);
1853       len         += len_si[proc];
1854     }
1855   }
1856 
1857   /* determine the number and length of messages to receive for coi and coj  */
1858   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1859   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1860 
1861   /* post the Irecv and Isend of coj */
1862   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1863   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1864   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1865   for (proc=0, k=0; proc<size; proc++) {
1866     if (!len_s[proc]) continue;
1867     i    = owners_co[proc];
1868     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1869     k++;
1870   }
1871 
1872   /* receives and sends of coj are complete */
1873   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1874   for (i=0; i<merge->nrecv; i++) {
1875     PetscMPIInt icompleted;
1876     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1877   }
1878   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1879   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1880 
1881   /* add received column indices into table to update Armax */
1882   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */
1883   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
1884     Jptr = buf_rj[k];
1885     for (j=0; j<merge->len_r[k]; j++) {
1886       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
1887     }
1888   }
1889   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
1890   /* 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); */
1891 
1892   /* send and recv coi */
1893   /*-------------------*/
1894   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1895   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1896   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1897   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1898   for (proc=0,k=0; proc<size; proc++) {
1899     if (!len_s[proc]) continue;
1900     /* form outgoing message for i-structure:
1901          buf_si[0]:                 nrows to be sent
1902                [1:nrows]:           row index (global)
1903                [nrows+1:2*nrows+1]: i-structure index
1904     */
1905     /*-------------------------------------------*/
1906     nrows       = len_si[proc]/2 - 1;
1907     buf_si_i    = buf_si + nrows+1;
1908     buf_si[0]   = nrows;
1909     buf_si_i[0] = 0;
1910     nrows       = 0;
1911     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1912       nzi               = coi[i+1] - coi[i];
1913       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1914       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1915       nrows++;
1916     }
1917     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1918     k++;
1919     buf_si += len_si[proc];
1920   }
1921   i = merge->nrecv;
1922   while (i--) {
1923     PetscMPIInt icompleted;
1924     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1925   }
1926   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1927   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1928   ierr = PetscFree(len_si);CHKERRQ(ierr);
1929   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1930   ierr = PetscFree(swaits);CHKERRQ(ierr);
1931   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1932   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1933 
1934   /* compute the local portion of C (mpi mat) */
1935   /*------------------------------------------*/
1936   /* allocate bi array and free space for accumulating nonzero column info */
1937   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1938   bi[0] = 0;
1939 
1940   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1941   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1942   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1943   current_space = free_space;
1944 
1945   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1946   for (k=0; k<merge->nrecv; k++) {
1947     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1948     nrows       = *buf_ri_k[k];
1949     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1950     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1951   }
1952 
1953   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1954   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1955   rmax = 0;
1956   for (i=0; i<pn; i++) {
1957     /* add pdt[i,:]*AP into lnk */
1958     pnz = pdti[i+1] - pdti[i];
1959     ptJ = pdtj + pdti[i];
1960     for (j=0; j<pnz; j++) {
1961       row  = ptJ[j];  /* row of AP == col of Pt */
1962       anz  = ai[row+1] - ai[row];
1963       Jptr = aj + ai[row];
1964       /* add non-zero cols of AP into the sorted linked list lnk */
1965       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1966     }
1967 
1968     /* add received col data into lnk */
1969     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1970       if (i == *nextrow[k]) { /* i-th row */
1971         nzi  = *(nextci[k]+1) - *nextci[k];
1972         Jptr = buf_rj[k] + *nextci[k];
1973         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1974         nextrow[k]++; nextci[k]++;
1975       }
1976     }
1977     nnz = lnk[0];
1978 
1979     /* if free space is not available, make more free space */
1980     if (current_space->local_remaining<nnz) {
1981       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1982       nspacedouble++;
1983     }
1984     /* copy data into free space, then initialize lnk */
1985     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1986     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1987 
1988     current_space->array           += nnz;
1989     current_space->local_used      += nnz;
1990     current_space->local_remaining -= nnz;
1991 
1992     bi[i+1] = bi[i] + nnz;
1993     if (nnz > rmax) rmax = nnz;
1994   }
1995   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1996 
1997   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1998   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1999   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
2000   if (afill_tmp > afill) afill = afill_tmp;
2001   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
2002   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
2003 
2004   ierr = MatDestroy(&POt);CHKERRQ(ierr);
2005   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
2006 
2007   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
2008   /*----------------------------------------------------------------------------------*/
2009   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
2010 
2011   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
2012   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
2013   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
2014   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
2015   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
2016   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
2017   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2018   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
2019   for (i=0; i<pn; i++) {
2020     row  = i + rstart;
2021     nnz  = bi[i+1] - bi[i];
2022     Jptr = bj + bi[i];
2023     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
2024   }
2025   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2026   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2027   ierr = PetscFree(vals);CHKERRQ(ierr);
2028 
2029   merge->bi        = bi;
2030   merge->bj        = bj;
2031   merge->coi       = coi;
2032   merge->coj       = coj;
2033   merge->buf_ri    = buf_ri;
2034   merge->buf_rj    = buf_rj;
2035   merge->owners_co = owners_co;
2036 
2037   /* attach the supporting struct to Cmpi for reuse */
2038   c = (Mat_MPIAIJ*)Cmpi->data;
2039 
2040   c->ptap     = ptap;
2041   ptap->api   = NULL;
2042   ptap->apj   = NULL;
2043   ptap->merge = merge;
2044   ptap->apa   = NULL;
2045   ptap->destroy   = Cmpi->ops->destroy;
2046   ptap->duplicate = Cmpi->ops->duplicate;
2047 
2048   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
2049   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
2050   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
2051 
2052   *C = Cmpi;
2053 #if defined(PETSC_USE_INFO)
2054   if (bi[pn] != 0) {
2055     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
2056     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
2057   } else {
2058     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
2059   }
2060 #endif
2061   PetscFunctionReturn(0);
2062 }
2063