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