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