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