xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 0db4d2e0165a0ea245cec0549c2de0bb7b39e2c0)
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 #include <petsc/private/sfimpl.h>
13 
14 #if defined(PETSC_HAVE_HYPRE)
15 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
16 #endif
17 
18 PETSC_INTERN PetscErrorCode MatProductSymbolic_AB_MPIAIJ_MPIAIJ(Mat C)
19 {
20   PetscErrorCode      ierr;
21   Mat_Product         *product = C->product;
22   Mat                 A=product->A,B=product->B;
23   MatProductAlgorithm alg=product->alg;
24   PetscReal           fill=product->fill;
25   PetscBool           flg;
26 
27   PetscFunctionBegin;
28   /* scalable */
29   ierr = PetscStrcmp(alg,"scalable",&flg);CHKERRQ(ierr);
30   if (flg) {
31     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
32     PetscFunctionReturn(0);
33   }
34 
35   /* nonscalable */
36   ierr = PetscStrcmp(alg,"nonscalable",&flg);CHKERRQ(ierr);
37   if (flg) {
38     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
39     PetscFunctionReturn(0);
40   }
41 
42   /* seqmpi */
43   ierr = PetscStrcmp(alg,"seqmpi",&flg);CHKERRQ(ierr);
44   if (flg) {
45     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(A,B,fill,C);CHKERRQ(ierr);
46     PetscFunctionReturn(0);
47   }
48 
49 #if defined(PETSC_HAVE_HYPRE)
50   ierr = PetscStrcmp(alg,"hypre",&flg);CHKERRQ(ierr);
51   if (flg) {
52     ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
53     PetscFunctionReturn(0);
54   }
55 #endif
56   SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Mat Product Algorithm is not supported");
57 }
58 
59 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(void *data)
60 {
61   PetscErrorCode ierr;
62   Mat_APMPI      *ptap = (Mat_APMPI*)data;
63 
64   PetscFunctionBegin;
65   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
66   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
67   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
68   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
69   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
70   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
71   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
72   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
73   ierr = PetscFree(ptap);CHKERRQ(ierr);
74   PetscFunctionReturn(0);
75 }
76 
77 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
78 {
79   PetscErrorCode ierr;
80   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
81   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
82   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
83   PetscScalar    *cda=cd->a,*coa=co->a;
84   Mat_SeqAIJ     *p_loc,*p_oth;
85   PetscScalar    *apa,*ca;
86   PetscInt       cm =C->rmap->n;
87   Mat_APMPI      *ptap;
88   PetscInt       *api,*apj,*apJ,i,k;
89   PetscInt       cstart=C->cmap->rstart;
90   PetscInt       cdnz,conz,k0,k1;
91   MPI_Comm       comm;
92   PetscMPIInt    size;
93 
94   PetscFunctionBegin;
95   MatCheckProduct(C,3);
96   ptap = (Mat_APMPI*)C->product->data;
97   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
98   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
99   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
100 
101   if (!ptap->P_oth && size>1) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatProductClear()");
102 
103   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
104   /*-----------------------------------------------------*/
105   /* update numerical values of P_oth and P_loc */
106   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
107   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
108 
109   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
110   /*----------------------------------------------------------*/
111   /* get data from symbolic products */
112   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
113   p_oth = NULL;
114   if (size >1) {
115     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
116   }
117 
118   /* get apa for storing dense row A[i,:]*P */
119   apa = ptap->apa;
120 
121   api = ptap->api;
122   apj = ptap->apj;
123   for (i=0; i<cm; i++) {
124     /* compute apa = A[i,:]*P */
125     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
126 
127     /* set values in C */
128     apJ  = apj + api[i];
129     cdnz = cd->i[i+1] - cd->i[i];
130     conz = co->i[i+1] - co->i[i];
131 
132     /* 1st off-diagonal part of C */
133     ca = coa + co->i[i];
134     k  = 0;
135     for (k0=0; k0<conz; k0++) {
136       if (apJ[k] >= cstart) break;
137       ca[k0]      = apa[apJ[k]];
138       apa[apJ[k++]] = 0.0;
139     }
140 
141     /* diagonal part of C */
142     ca = cda + cd->i[i];
143     for (k1=0; k1<cdnz; k1++) {
144       ca[k1]      = apa[apJ[k]];
145       apa[apJ[k++]] = 0.0;
146     }
147 
148     /* 2nd off-diagonal part of C */
149     ca = coa + co->i[i];
150     for (; k0<conz; k0++) {
151       ca[k0]      = apa[apJ[k]];
152       apa[apJ[k++]] = 0.0;
153     }
154   }
155   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157   PetscFunctionReturn(0);
158 }
159 
160 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat C)
161 {
162   PetscErrorCode     ierr;
163   MPI_Comm           comm;
164   PetscMPIInt        size;
165   Mat_APMPI          *ptap;
166   PetscFreeSpaceList free_space=NULL,current_space=NULL;
167   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
168   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
169   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
170   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
171   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
172   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
173   PetscBT            lnkbt;
174   PetscReal          afill;
175   MatType            mtype;
176 
177   PetscFunctionBegin;
178   MatCheckProduct(C,4);
179   if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
180   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
181   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
182 
183   /* create struct Mat_APMPI and attached it to C later */
184   ierr = PetscNew(&ptap);CHKERRQ(ierr);
185 
186   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
187   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
188 
189   /* get P_loc by taking all local rows of P */
190   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
191 
192   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
193   pi_loc = p_loc->i; pj_loc = p_loc->j;
194   if (size > 1) {
195     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
196     pi_oth = p_oth->i; pj_oth = p_oth->j;
197   } else {
198     p_oth = NULL;
199     pi_oth = NULL; pj_oth = NULL;
200   }
201 
202   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
203   /*-------------------------------------------------------------------*/
204   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
205   ptap->api = api;
206   api[0]    = 0;
207 
208   /* create and initialize a linked list */
209   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
210 
211   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
212   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
213   current_space = free_space;
214 
215   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
216   for (i=0; i<am; i++) {
217     /* diagonal portion of A */
218     nzi = adi[i+1] - adi[i];
219     for (j=0; j<nzi; j++) {
220       row  = *adj++;
221       pnz  = pi_loc[row+1] - pi_loc[row];
222       Jptr = pj_loc + pi_loc[row];
223       /* add non-zero cols of P into the sorted linked list lnk */
224       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
225     }
226     /* off-diagonal portion of A */
227     nzi = aoi[i+1] - aoi[i];
228     for (j=0; j<nzi; j++) {
229       row  = *aoj++;
230       pnz  = pi_oth[row+1] - pi_oth[row];
231       Jptr = pj_oth + pi_oth[row];
232       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
233     }
234     /* add possible missing diagonal entry */
235     if (C->force_diagonals) {
236       j = i + rstart; /* column index */
237       ierr = PetscLLCondensedAddSorted(1,&j,lnk,lnkbt);CHKERRQ(ierr);
238     }
239 
240     apnz     = lnk[0];
241     api[i+1] = api[i] + apnz;
242 
243     /* if free space is not available, double the total space in the list */
244     if (current_space->local_remaining<apnz) {
245       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
246       nspacedouble++;
247     }
248 
249     /* Copy data into free space, then initialize lnk */
250     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
251     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
252 
253     current_space->array           += apnz;
254     current_space->local_used      += apnz;
255     current_space->local_remaining -= apnz;
256   }
257 
258   /* Allocate space for apj, initialize apj, and */
259   /* destroy list of free space and other temporary array(s) */
260   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
261   apj  = ptap->apj;
262   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
263   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
264 
265   /* malloc apa to store dense row A[i,:]*P */
266   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
267 
268   /* set and assemble symbolic parallel matrix C */
269   /*---------------------------------------------*/
270   ierr = MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
271   ierr = MatSetBlockSizesFromMats(C,A,P);CHKERRQ(ierr);
272 
273   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
274   ierr = MatSetType(C,mtype);CHKERRQ(ierr);
275   ierr = MatMPIAIJSetPreallocation(C,0,dnz,0,onz);CHKERRQ(ierr);
276   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
277 
278   ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);CHKERRQ(ierr);
279   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
280   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
281   ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
282 
283   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
284   C->ops->productnumeric = MatProductNumeric_AB;
285 
286   /* attach the supporting struct to C for reuse */
287   C->product->data    = ptap;
288   C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;
289 
290   /* set MatInfo */
291   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
292   if (afill < 1.0) afill = 1.0;
293   C->info.mallocs           = nspacedouble;
294   C->info.fill_ratio_given  = fill;
295   C->info.fill_ratio_needed = afill;
296 
297 #if defined(PETSC_USE_INFO)
298   if (api[am]) {
299     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
300     ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
301   } else {
302     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
303   }
304 #endif
305   PetscFunctionReturn(0);
306 }
307 
308 /* ------------------------------------------------------- */
309 static PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat);
310 static PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
311 
312 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AB(Mat C)
313 {
314   Mat_Product *product = C->product;
315   Mat         A = product->A,B=product->B;
316 
317   PetscFunctionBegin;
318   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend)
319     SETERRQ4(PETSC_COMM_SELF,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);
320 
321   C->ops->matmultsymbolic = MatMatMultSymbolic_MPIAIJ_MPIDense;
322   C->ops->productsymbolic = MatProductSymbolic_AB;
323   PetscFunctionReturn(0);
324 }
325 /* -------------------------------------------------------------------- */
326 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(Mat C)
327 {
328   Mat_Product *product = C->product;
329   Mat         A = product->A,B=product->B;
330 
331   PetscFunctionBegin;
332   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend)
333     SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);
334 
335   C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_MPIAIJ_MPIDense;
336   C->ops->productsymbolic          = MatProductSymbolic_AtB;
337   PetscFunctionReturn(0);
338 }
339 
340 /* --------------------------------------------------------------------- */
341 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat C)
342 {
343   PetscErrorCode ierr;
344   Mat_Product    *product = C->product;
345 
346   PetscFunctionBegin;
347   switch (product->type) {
348   case MATPRODUCT_AB:
349     ierr = MatProductSetFromOptions_MPIAIJ_MPIDense_AB(C);CHKERRQ(ierr);
350     break;
351   case MATPRODUCT_AtB:
352     ierr = MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(C);CHKERRQ(ierr);
353     break;
354   default:
355     break;
356   }
357   PetscFunctionReturn(0);
358 }
359 /* ------------------------------------------------------- */
360 
361 typedef struct {
362   Mat          workB,workB1;
363   MPI_Request  *rwaits,*swaits;
364   PetscInt     nsends,nrecvs;
365   MPI_Datatype *stype,*rtype;
366   PetscInt     blda;
367 } MPIAIJ_MPIDense;
368 
369 PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
370 {
371   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*)ctx;
372   PetscErrorCode  ierr;
373   PetscInt        i;
374 
375   PetscFunctionBegin;
376   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
377   ierr = MatDestroy(&contents->workB1);CHKERRQ(ierr);
378   for (i=0; i<contents->nsends; i++) {
379     ierr = MPI_Type_free(&contents->stype[i]);CHKERRMPI(ierr);
380   }
381   for (i=0; i<contents->nrecvs; i++) {
382     ierr = MPI_Type_free(&contents->rtype[i]);CHKERRMPI(ierr);
383   }
384   ierr = PetscFree4(contents->stype,contents->rtype,contents->rwaits,contents->swaits);CHKERRQ(ierr);
385   ierr = PetscFree(contents);CHKERRQ(ierr);
386   PetscFunctionReturn(0);
387 }
388 
389 static PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat C)
390 {
391   PetscErrorCode  ierr;
392   Mat_MPIAIJ      *aij=(Mat_MPIAIJ*)A->data;
393   PetscInt        nz=aij->B->cmap->n,nsends,nrecvs,i,nrows_to,j,blda,clda;
394   MPIAIJ_MPIDense *contents;
395   VecScatter      ctx=aij->Mvctx;
396   PetscInt        Am=A->rmap->n,Bm=B->rmap->n,BN=B->cmap->N,Bbn,Bbn1,bs,nrows_from,numBb;
397   MPI_Comm        comm;
398   MPI_Datatype    type1,*stype,*rtype;
399   const PetscInt  *sindices,*sstarts,*rstarts;
400   PetscMPIInt     *disp;
401   PetscBool       cisdense;
402 
403   PetscFunctionBegin;
404   MatCheckProduct(C,4);
405   if (C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
406   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
407   ierr = PetscObjectBaseTypeCompare((PetscObject)C,MATMPIDENSE,&cisdense);CHKERRQ(ierr);
408   if (!cisdense) {
409     ierr = MatSetType(C,((PetscObject)B)->type_name);CHKERRQ(ierr);
410   }
411   ierr = MatSetSizes(C,Am,B->cmap->n,A->rmap->N,BN);CHKERRQ(ierr);
412   ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr);
413   ierr = MatSetUp(C);CHKERRQ(ierr);
414   ierr = MatDenseGetLDA(B,&blda);CHKERRQ(ierr);
415   ierr = MatDenseGetLDA(C,&clda);CHKERRQ(ierr);
416   ierr = PetscNew(&contents);CHKERRQ(ierr);
417 
418   ierr = VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,NULL,NULL);CHKERRQ(ierr);
419   ierr = VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,NULL,NULL);CHKERRQ(ierr);
420 
421   /* Create column block of B and C for memory scalability when BN is too large */
422   /* Estimate Bbn, column size of Bb */
423   if (nz) {
424     Bbn1 = 2*Am*BN/nz;
425   } else Bbn1 = BN;
426 
427   bs = PetscAbs(B->cmap->bs);
428   Bbn1 = Bbn1/bs *bs; /* Bbn1 is a multiple of bs */
429   if (Bbn1 > BN) Bbn1 = BN;
430   ierr = MPI_Allreduce(&Bbn1,&Bbn,1,MPIU_INT,MPI_MAX,comm);CHKERRMPI(ierr);
431 
432   /* Enable runtime option for Bbn */
433   ierr = PetscOptionsBegin(comm,((PetscObject)C)->prefix,"MatMatMult","Mat");CHKERRQ(ierr);
434   ierr = PetscOptionsInt("-matmatmult_Bbn","Number of columns in Bb","MatMatMult",Bbn,&Bbn,NULL);CHKERRQ(ierr);
435   ierr = PetscOptionsEnd();CHKERRQ(ierr);
436   Bbn  = PetscMin(Bbn,BN);
437 
438   if (Bbn > 0 && Bbn < BN) {
439     numBb = BN/Bbn;
440     Bbn1 = BN - numBb*Bbn;
441   } else numBb = 0;
442 
443   if (numBb) {
444     ierr = PetscInfo3(C,"use Bb, BN=%D, Bbn=%D; numBb=%D\n",BN,Bbn,numBb);CHKERRQ(ierr);
445     if (Bbn1) { /* Create workB1 for the remaining columns */
446       ierr = PetscInfo2(C,"use Bb1, BN=%D, Bbn1=%D\n",BN,Bbn1);CHKERRQ(ierr);
447       /* Create work matrix used to store off processor rows of B needed for local product */
448       ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,Bbn1,NULL,&contents->workB1);CHKERRQ(ierr);
449     } else contents->workB1 = NULL;
450   }
451 
452   /* Create work matrix used to store off processor rows of B needed for local product */
453   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,Bbn,NULL,&contents->workB);CHKERRQ(ierr);
454 
455   /* Use MPI derived data type to reduce memory required by the send/recv buffers */
456   ierr = PetscMalloc4(nsends,&stype,nrecvs,&rtype,nrecvs,&contents->rwaits,nsends,&contents->swaits);CHKERRQ(ierr);
457   contents->stype  = stype;
458   contents->nsends = nsends;
459 
460   contents->rtype  = rtype;
461   contents->nrecvs = nrecvs;
462   contents->blda   = blda;
463 
464   ierr = PetscMalloc1(Bm+1,&disp);CHKERRQ(ierr);
465   for (i=0; i<nsends; i++) {
466     nrows_to = sstarts[i+1]-sstarts[i];
467     for (j=0; j<nrows_to; j++){
468       disp[j] = sindices[sstarts[i]+j]; /* rowB to be sent */
469     }
470     ierr = MPI_Type_create_indexed_block(nrows_to,1,(const PetscMPIInt *)disp,MPIU_SCALAR,&type1);CHKERRMPI(ierr);
471 
472     ierr = MPI_Type_create_resized(type1,0,blda*sizeof(PetscScalar),&stype[i]);CHKERRMPI(ierr);
473     ierr = MPI_Type_commit(&stype[i]);CHKERRMPI(ierr);
474     ierr = MPI_Type_free(&type1);CHKERRMPI(ierr);
475   }
476 
477   for (i=0; i<nrecvs; i++) {
478     /* received values from a process form a (nrows_from x Bbn) row block in workB (column-wise) */
479     nrows_from = rstarts[i+1]-rstarts[i];
480     disp[0] = 0;
481     ierr = MPI_Type_create_indexed_block(1, nrows_from, (const PetscMPIInt *)disp, MPIU_SCALAR, &type1);CHKERRMPI(ierr);
482     ierr = MPI_Type_create_resized(type1, 0, nz*sizeof(PetscScalar), &rtype[i]);CHKERRMPI(ierr);
483     ierr = MPI_Type_commit(&rtype[i]);CHKERRMPI(ierr);
484     ierr = MPI_Type_free(&type1);CHKERRMPI(ierr);
485   }
486 
487   ierr = PetscFree(disp);CHKERRQ(ierr);
488   ierr = VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,NULL,NULL);CHKERRQ(ierr);
489   ierr = VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,NULL,NULL);CHKERRQ(ierr);
490   ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
491   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
492   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
493   ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
494 
495   C->product->data = contents;
496   C->product->destroy = MatMPIAIJ_MPIDenseDestroy;
497   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
498   PetscFunctionReturn(0);
499 }
500 
501 PETSC_INTERN PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat,const PetscBool);
502 /*
503     Performs an efficient scatter on the rows of B needed by this process; this is
504     a modification of the VecScatterBegin_() routines.
505 
506     Input: Bbidx = 0: B = Bb
507                  = 1: B = Bb1, see MatMatMultSymbolic_MPIAIJ_MPIDense()
508 */
509 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,PetscInt Bbidx,Mat C,Mat *outworkB)
510 {
511   Mat_MPIAIJ        *aij = (Mat_MPIAIJ*)A->data;
512   PetscErrorCode    ierr;
513   const PetscScalar *b;
514   PetscScalar       *rvalues;
515   VecScatter        ctx = aij->Mvctx;
516   const PetscInt    *sindices,*sstarts,*rstarts;
517   const PetscMPIInt *sprocs,*rprocs;
518   PetscInt          i,nsends,nrecvs;
519   MPI_Request       *swaits,*rwaits;
520   MPI_Comm          comm;
521   PetscMPIInt       tag=((PetscObject)ctx)->tag,ncols=B->cmap->N,nrows=aij->B->cmap->n,nsends_mpi,nrecvs_mpi;
522   MPIAIJ_MPIDense   *contents;
523   Mat               workB;
524   MPI_Datatype      *stype,*rtype;
525   PetscInt          blda;
526 
527   PetscFunctionBegin;
528   MatCheckProduct(C,4);
529   if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
530   contents = (MPIAIJ_MPIDense*)C->product->data;
531   ierr = VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL/*bs*/);CHKERRQ(ierr);
532   ierr = VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL/*bs*/);CHKERRQ(ierr);
533   ierr = PetscMPIIntCast(nsends,&nsends_mpi);CHKERRQ(ierr);
534   ierr = PetscMPIIntCast(nrecvs,&nrecvs_mpi);CHKERRQ(ierr);
535   if (Bbidx == 0) {
536     workB = *outworkB = contents->workB;
537   } else {
538     workB = *outworkB = contents->workB1;
539   }
540   if (nrows != workB->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",workB->cmap->n,nrows);
541   swaits = contents->swaits;
542   rwaits = contents->rwaits;
543 
544   ierr = MatDenseGetArrayRead(B,&b);CHKERRQ(ierr);
545   ierr = MatDenseGetLDA(B,&blda);CHKERRQ(ierr);
546   if (blda != contents->blda) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot reuse an input matrix with lda %D != %D",blda,contents->blda);
547   ierr = MatDenseGetArray(workB,&rvalues);CHKERRQ(ierr);
548 
549   /* Post recv, use MPI derived data type to save memory */
550   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
551   rtype = contents->rtype;
552   for (i=0; i<nrecvs; i++) {
553     ierr = MPI_Irecv(rvalues+(rstarts[i]-rstarts[0]),ncols,rtype[i],rprocs[i],tag,comm,rwaits+i);CHKERRMPI(ierr);
554   }
555 
556   stype = contents->stype;
557   for (i=0; i<nsends; i++) {
558     ierr = MPI_Isend(b,ncols,stype[i],sprocs[i],tag,comm,swaits+i);CHKERRMPI(ierr);
559   }
560 
561   if (nrecvs) {ierr = MPI_Waitall(nrecvs_mpi,rwaits,MPI_STATUSES_IGNORE);CHKERRMPI(ierr);}
562   if (nsends) {ierr = MPI_Waitall(nsends_mpi,swaits,MPI_STATUSES_IGNORE);CHKERRMPI(ierr);}
563 
564   ierr = VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL);CHKERRQ(ierr);
565   ierr = VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL);CHKERRQ(ierr);
566   ierr = MatDenseRestoreArrayRead(B,&b);CHKERRQ(ierr);
567   ierr = MatDenseRestoreArray(workB,&rvalues);CHKERRQ(ierr);
568   PetscFunctionReturn(0);
569 }
570 
571 static PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
572 {
573   PetscErrorCode  ierr;
574   Mat_MPIAIJ      *aij    = (Mat_MPIAIJ*)A->data;
575   Mat_MPIDense    *bdense = (Mat_MPIDense*)B->data;
576   Mat_MPIDense    *cdense = (Mat_MPIDense*)C->data;
577   Mat             workB;
578   MPIAIJ_MPIDense *contents;
579 
580   PetscFunctionBegin;
581   MatCheckProduct(C,3);
582   if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
583   contents = (MPIAIJ_MPIDense*)C->product->data;
584   /* diagonal block of A times all local rows of B */
585   /* TODO: this calls a symbolic multiplication every time, which could be avoided */
586   ierr = MatMatMult(aij->A,bdense->A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&cdense->A);CHKERRQ(ierr);
587   if (contents->workB->cmap->n == B->cmap->N) {
588     /* get off processor parts of B needed to complete C=A*B */
589     ierr = MatMPIDenseScatter(A,B,0,C,&workB);CHKERRQ(ierr);
590 
591     /* off-diagonal block of A times nonlocal rows of B */
592     ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A,PETSC_TRUE);CHKERRQ(ierr);
593   } else {
594     Mat      Bb,Cb;
595     PetscInt BN=B->cmap->N,n=contents->workB->cmap->n,i;
596 
597     for (i=0; i<BN; i+=n) {
598       ierr = MatDenseGetSubMatrix(B,i,PetscMin(i+n,BN),&Bb);CHKERRQ(ierr);
599       ierr = MatDenseGetSubMatrix(C,i,PetscMin(i+n,BN),&Cb);CHKERRQ(ierr);
600 
601       /* get off processor parts of B needed to complete C=A*B */
602       ierr = MatMPIDenseScatter(A,Bb,i+n>BN,C,&workB);CHKERRQ(ierr);
603 
604       /* off-diagonal block of A times nonlocal rows of B */
605       cdense = (Mat_MPIDense*)Cb->data;
606       ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A,PETSC_TRUE);CHKERRQ(ierr);
607 
608       ierr = MatDenseRestoreSubMatrix(B,&Bb);CHKERRQ(ierr);
609       ierr = MatDenseRestoreSubMatrix(C,&Cb);CHKERRQ(ierr);
610     }
611   }
612   PetscFunctionReturn(0);
613 }
614 
615 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
616 {
617   PetscErrorCode ierr;
618   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
619   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
620   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
621   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
622   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
623   Mat_SeqAIJ     *p_loc,*p_oth;
624   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
625   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
626   PetscInt       cm    = C->rmap->n,anz,pnz;
627   Mat_APMPI      *ptap;
628   PetscScalar    *apa_sparse;
629   PetscInt       *api,*apj,*apJ,i,j,k,row;
630   PetscInt       cstart = C->cmap->rstart;
631   PetscInt       cdnz,conz,k0,k1,nextp;
632   MPI_Comm       comm;
633   PetscMPIInt    size;
634 
635   PetscFunctionBegin;
636   MatCheckProduct(C,3);
637   ptap = (Mat_APMPI*)C->product->data;
638   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
639   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
640   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
641   if (!ptap->P_oth && size>1) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatProductClear()");
642 
643   apa_sparse = ptap->apa;
644 
645   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
646   /*-----------------------------------------------------*/
647   /* update numerical values of P_oth and P_loc */
648   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
649   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
650 
651   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
652   /*----------------------------------------------------------*/
653   /* get data from symbolic products */
654   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
655   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
656   if (size >1) {
657     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
658     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
659   } else {
660     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
661   }
662 
663   api = ptap->api;
664   apj = ptap->apj;
665   for (i=0; i<cm; i++) {
666     apJ = apj + api[i];
667 
668     /* diagonal portion of A */
669     anz = adi[i+1] - adi[i];
670     adj = ad->j + adi[i];
671     ada = ad->a + adi[i];
672     for (j=0; j<anz; j++) {
673       row = adj[j];
674       pnz = pi_loc[row+1] - pi_loc[row];
675       pj  = pj_loc + pi_loc[row];
676       pa  = pa_loc + pi_loc[row];
677       /* perform sparse axpy */
678       valtmp = ada[j];
679       nextp  = 0;
680       for (k=0; nextp<pnz; k++) {
681         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
682           apa_sparse[k] += valtmp*pa[nextp++];
683         }
684       }
685       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
686     }
687 
688     /* off-diagonal portion of A */
689     anz = aoi[i+1] - aoi[i];
690     aoj = ao->j + aoi[i];
691     aoa = ao->a + aoi[i];
692     for (j=0; j<anz; j++) {
693       row = aoj[j];
694       pnz = pi_oth[row+1] - pi_oth[row];
695       pj  = pj_oth + pi_oth[row];
696       pa  = pa_oth + pi_oth[row];
697       /* perform sparse axpy */
698       valtmp = aoa[j];
699       nextp  = 0;
700       for (k=0; nextp<pnz; k++) {
701         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
702           apa_sparse[k] += valtmp*pa[nextp++];
703         }
704       }
705       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
706     }
707 
708     /* set values in C */
709     cdnz = cd->i[i+1] - cd->i[i];
710     conz = co->i[i+1] - co->i[i];
711 
712     /* 1st off-diagonal part of C */
713     ca = coa + co->i[i];
714     k  = 0;
715     for (k0=0; k0<conz; k0++) {
716       if (apJ[k] >= cstart) break;
717       ca[k0]        = apa_sparse[k];
718       apa_sparse[k] = 0.0;
719       k++;
720     }
721 
722     /* diagonal part of C */
723     ca = cda + cd->i[i];
724     for (k1=0; k1<cdnz; k1++) {
725       ca[k1]        = apa_sparse[k];
726       apa_sparse[k] = 0.0;
727       k++;
728     }
729 
730     /* 2nd off-diagonal part of C */
731     ca = coa + co->i[i];
732     for (; k0<conz; k0++) {
733       ca[k0]        = apa_sparse[k];
734       apa_sparse[k] = 0.0;
735       k++;
736     }
737   }
738   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
739   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
740   PetscFunctionReturn(0);
741 }
742 
743 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
744 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat C)
745 {
746   PetscErrorCode     ierr;
747   MPI_Comm           comm;
748   PetscMPIInt        size;
749   Mat_APMPI          *ptap;
750   PetscFreeSpaceList free_space = NULL,current_space=NULL;
751   Mat_MPIAIJ         *a  = (Mat_MPIAIJ*)A->data;
752   Mat_SeqAIJ         *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
753   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
754   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
755   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=1;
756   PetscInt           am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20;
757   PetscReal          afill;
758   MatType            mtype;
759 
760   PetscFunctionBegin;
761   MatCheckProduct(C,4);
762   if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
763   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
764   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
765 
766   /* create struct Mat_APMPI and attached it to C later */
767   ierr = PetscNew(&ptap);CHKERRQ(ierr);
768 
769   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
770   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
771 
772   /* get P_loc by taking all local rows of P */
773   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
774 
775   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
776   pi_loc = p_loc->i; pj_loc = p_loc->j;
777   if (size > 1) {
778     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
779     pi_oth = p_oth->i; pj_oth = p_oth->j;
780   } else {
781     p_oth  = NULL;
782     pi_oth = NULL; pj_oth = NULL;
783   }
784 
785   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
786   /*-------------------------------------------------------------------*/
787   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
788   ptap->api = api;
789   api[0]    = 0;
790 
791   ierr = PetscLLCondensedCreate_Scalable(lsize,&lnk);CHKERRQ(ierr);
792 
793   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
794   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
795   current_space = free_space;
796   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
797   for (i=0; i<am; i++) {
798     /* diagonal portion of A */
799     nzi = adi[i+1] - adi[i];
800     for (j=0; j<nzi; j++) {
801       row  = *adj++;
802       pnz  = pi_loc[row+1] - pi_loc[row];
803       Jptr = pj_loc + pi_loc[row];
804       /* Expand list if it is not long enough */
805       if (pnz+apnz_max > lsize) {
806         lsize = pnz+apnz_max;
807         ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr);
808       }
809       /* add non-zero cols of P into the sorted linked list lnk */
810       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
811       apnz     = *lnk; /* The first element in the list is the number of items in the list */
812       api[i+1] = api[i] + apnz;
813       if (apnz > apnz_max) apnz_max = apnz + 1; /* '1' for diagonal entry */
814     }
815     /* off-diagonal portion of A */
816     nzi = aoi[i+1] - aoi[i];
817     for (j=0; j<nzi; j++) {
818       row  = *aoj++;
819       pnz  = pi_oth[row+1] - pi_oth[row];
820       Jptr = pj_oth + pi_oth[row];
821       /* Expand list if it is not long enough */
822       if (pnz+apnz_max > lsize) {
823         lsize = pnz + apnz_max;
824         ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr);
825       }
826       /* add non-zero cols of P into the sorted linked list lnk */
827       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
828       apnz     = *lnk;  /* The first element in the list is the number of items in the list */
829       api[i+1] = api[i] + apnz;
830       if (apnz > apnz_max) apnz_max = apnz + 1; /* '1' for diagonal entry */
831     }
832 
833     /* add missing diagonal entry */
834     if (C->force_diagonals) {
835       j = i + rstart; /* column index */
836       ierr = PetscLLCondensedAddSorted_Scalable(1,&j,lnk);CHKERRQ(ierr);
837     }
838 
839     apnz     = *lnk;
840     api[i+1] = api[i] + apnz;
841     if (apnz > apnz_max) apnz_max = apnz;
842 
843     /* if free space is not available, double the total space in the list */
844     if (current_space->local_remaining<apnz) {
845       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
846       nspacedouble++;
847     }
848 
849     /* Copy data into free space, then initialize lnk */
850     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
851     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
852 
853     current_space->array           += apnz;
854     current_space->local_used      += apnz;
855     current_space->local_remaining -= apnz;
856   }
857 
858   /* Allocate space for apj, initialize apj, and */
859   /* destroy list of free space and other temporary array(s) */
860   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
861   apj  = ptap->apj;
862   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
863   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
864 
865   /* create and assemble symbolic parallel matrix C */
866   /*----------------------------------------------------*/
867   ierr = MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
868   ierr = MatSetBlockSizesFromMats(C,A,P);CHKERRQ(ierr);
869   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
870   ierr = MatSetType(C,mtype);CHKERRQ(ierr);
871   ierr = MatMPIAIJSetPreallocation(C,0,dnz,0,onz);CHKERRQ(ierr);
872   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
873 
874   /* malloc apa for assembly C */
875   ierr = PetscCalloc1(apnz_max,&ptap->apa);CHKERRQ(ierr);
876 
877   ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);CHKERRQ(ierr);
878   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
879   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
880   ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
881 
882   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
883   C->ops->productnumeric = MatProductNumeric_AB;
884 
885   /* attach the supporting struct to C for reuse */
886   C->product->data    = ptap;
887   C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;
888 
889   /* set MatInfo */
890   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
891   if (afill < 1.0) afill = 1.0;
892   C->info.mallocs           = nspacedouble;
893   C->info.fill_ratio_given  = fill;
894   C->info.fill_ratio_needed = afill;
895 
896 #if defined(PETSC_USE_INFO)
897   if (api[am]) {
898     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
899     ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
900   } else {
901     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
902   }
903 #endif
904   PetscFunctionReturn(0);
905 }
906 
907 /* This function is needed for the seqMPI matrix-matrix multiplication.  */
908 /* Three input arrays are merged to one output array. The size of the    */
909 /* output array is also output. Duplicate entries only show up once.     */
910 static void Merge3SortedArrays(PetscInt  size1, PetscInt *in1,
911                                PetscInt  size2, PetscInt *in2,
912                                PetscInt  size3, PetscInt *in3,
913                                PetscInt *size4, PetscInt *out)
914 {
915   int i = 0, j = 0, k = 0, l = 0;
916 
917   /* Traverse all three arrays */
918   while (i<size1 && j<size2 && k<size3) {
919     if (in1[i] < in2[j] && in1[i] < in3[k]) {
920       out[l++] = in1[i++];
921     }
922     else if (in2[j] < in1[i] && in2[j] < in3[k]) {
923       out[l++] = in2[j++];
924     }
925     else if (in3[k] < in1[i] && in3[k] < in2[j]) {
926       out[l++] = in3[k++];
927     }
928     else if (in1[i] == in2[j] && in1[i] < in3[k]) {
929       out[l++] = in1[i];
930       i++, j++;
931     }
932     else if (in1[i] == in3[k] && in1[i] < in2[j]) {
933       out[l++] = in1[i];
934       i++, k++;
935     }
936     else if (in3[k] == in2[j] && in2[j] < in1[i])  {
937       out[l++] = in2[j];
938       k++, j++;
939     }
940     else if (in1[i] == in2[j] && in1[i] == in3[k]) {
941       out[l++] = in1[i];
942       i++, j++, k++;
943     }
944   }
945 
946   /* Traverse two remaining arrays */
947   while (i<size1 && j<size2) {
948     if (in1[i] < in2[j]) {
949       out[l++] = in1[i++];
950     }
951     else if (in1[i] > in2[j]) {
952       out[l++] = in2[j++];
953     }
954     else {
955       out[l++] = in1[i];
956       i++, j++;
957     }
958   }
959 
960   while (i<size1 && k<size3) {
961     if (in1[i] < in3[k]) {
962       out[l++] = in1[i++];
963     }
964     else if (in1[i] > in3[k]) {
965       out[l++] = in3[k++];
966     }
967     else {
968       out[l++] = in1[i];
969       i++, k++;
970     }
971   }
972 
973   while (k<size3 && j<size2)  {
974     if (in3[k] < in2[j]) {
975       out[l++] = in3[k++];
976     }
977     else if (in3[k] > in2[j]) {
978       out[l++] = in2[j++];
979     }
980     else {
981       out[l++] = in3[k];
982       k++, j++;
983     }
984   }
985 
986   /* Traverse one remaining array */
987   while (i<size1) out[l++] = in1[i++];
988   while (j<size2) out[l++] = in2[j++];
989   while (k<size3) out[l++] = in3[k++];
990 
991   *size4 = l;
992 }
993 
994 /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and  */
995 /* adds up the products. Two of these three multiplications are performed with existing (sequential)      */
996 /* matrix-matrix multiplications.  */
997 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat C)
998 {
999   PetscErrorCode     ierr;
1000   MPI_Comm           comm;
1001   PetscMPIInt        size;
1002   Mat_APMPI          *ptap;
1003   PetscFreeSpaceList free_space_diag=NULL, current_space=NULL;
1004   Mat_MPIAIJ         *a  =(Mat_MPIAIJ*)A->data;
1005   Mat_SeqAIJ         *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc;
1006   Mat_MPIAIJ         *p  =(Mat_MPIAIJ*)P->data;
1007   Mat_SeqAIJ         *adpd_seq, *p_off, *aopoth_seq;
1008   PetscInt           adponz, adpdnz;
1009   PetscInt           *pi_loc,*dnz,*onz;
1010   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart;
1011   PetscInt           *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi,
1012                      *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj;
1013   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend;
1014   PetscBT            lnkbt;
1015   PetscReal          afill;
1016   PetscMPIInt        rank;
1017   Mat                adpd, aopoth;
1018   MatType            mtype;
1019   const char         *prefix;
1020 
1021   PetscFunctionBegin;
1022   MatCheckProduct(C,4);
1023   if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
1024   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1025   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
1026   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
1027   ierr = MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend);CHKERRQ(ierr);
1028 
1029   /* create struct Mat_APMPI and attached it to C later */
1030   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1031 
1032   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1033   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1034 
1035   /* get P_loc by taking all local rows of P */
1036   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1037 
1038 
1039   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1040   pi_loc = p_loc->i;
1041 
1042   /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */
1043   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
1044   ierr      = PetscMalloc1(am+2,&adpoi);CHKERRQ(ierr);
1045 
1046   adpoi[0]    = 0;
1047   ptap->api = api;
1048   api[0] = 0;
1049 
1050   /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */
1051   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
1052   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
1053 
1054   /* Symbolic calc of A_loc_diag * P_loc_diag */
1055   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1056   ierr = MatProductCreate(a->A,p->A,NULL,&adpd);CHKERRQ(ierr);
1057   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1058   ierr = MatSetOptionsPrefix(adpd,prefix);CHKERRQ(ierr);
1059   ierr = MatAppendOptionsPrefix(adpd,"inner_diag_");CHKERRQ(ierr);
1060 
1061   ierr = MatProductSetType(adpd,MATPRODUCT_AB);CHKERRQ(ierr);
1062   ierr = MatProductSetAlgorithm(adpd,"sorted");CHKERRQ(ierr);
1063   ierr = MatProductSetFill(adpd,fill);CHKERRQ(ierr);
1064   ierr = MatProductSetFromOptions(adpd);CHKERRQ(ierr);
1065 
1066   adpd->force_diagonals = C->force_diagonals;
1067   ierr = MatProductSymbolic(adpd);CHKERRQ(ierr);
1068 
1069   adpd_seq = (Mat_SeqAIJ*)((adpd)->data);
1070   adpdi = adpd_seq->i; adpdj = adpd_seq->j;
1071   p_off = (Mat_SeqAIJ*)((p->B)->data);
1072   poff_i = p_off->i; poff_j = p_off->j;
1073 
1074   /* j_temp stores indices of a result row before they are added to the linked list */
1075   ierr = PetscMalloc1(pN+2,&j_temp);CHKERRQ(ierr);
1076 
1077 
1078   /* Symbolic calc of the A_diag * p_loc_off */
1079   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
1080   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);CHKERRQ(ierr);
1081   current_space = free_space_diag;
1082 
1083   for (i=0; i<am; i++) {
1084     /* A_diag * P_loc_off */
1085     nzi = adi[i+1] - adi[i];
1086     for (j=0; j<nzi; j++) {
1087       row  = *adj++;
1088       pnz  = poff_i[row+1] - poff_i[row];
1089       Jptr = poff_j + poff_i[row];
1090       for (i1 = 0; i1 < pnz; i1++) {
1091         j_temp[i1] = p->garray[Jptr[i1]];
1092       }
1093       /* add non-zero cols of P into the sorted linked list lnk */
1094       ierr = PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);CHKERRQ(ierr);
1095     }
1096 
1097     adponz     = lnk[0];
1098     adpoi[i+1] = adpoi[i] + adponz;
1099 
1100     /* if free space is not available, double the total space in the list */
1101     if (current_space->local_remaining<adponz) {
1102       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1103       nspacedouble++;
1104     }
1105 
1106     /* Copy data into free space, then initialize lnk */
1107     ierr = PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1108 
1109     current_space->array           += adponz;
1110     current_space->local_used      += adponz;
1111     current_space->local_remaining -= adponz;
1112   }
1113 
1114   /* Symbolic calc of A_off * P_oth */
1115   ierr = MatSetOptionsPrefix(a->B,prefix);CHKERRQ(ierr);
1116   ierr = MatAppendOptionsPrefix(a->B,"inner_offdiag_");CHKERRQ(ierr);
1117   ierr = MatCreate(PETSC_COMM_SELF,&aopoth);CHKERRQ(ierr);
1118   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, aopoth);CHKERRQ(ierr);
1119   aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data);
1120   aopothi = aopoth_seq->i; aopothj = aopoth_seq->j;
1121 
1122   /* Allocate space for apj, adpj, aopj, ... */
1123   /* destroy lists of free space and other temporary array(s) */
1124 
1125   ierr = PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);CHKERRQ(ierr);
1126   ierr = PetscMalloc1(adpoi[am]+2, &adpoj);CHKERRQ(ierr);
1127 
1128   /* Copy from linked list to j-array */
1129   ierr = PetscFreeSpaceContiguous(&free_space_diag,adpoj);CHKERRQ(ierr);
1130   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1131 
1132   adpoJ = adpoj;
1133   adpdJ = adpdj;
1134   aopJ = aopothj;
1135   apj  = ptap->apj;
1136   apJ = apj; /* still empty */
1137 
1138   /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1139   /* A_diag * P_loc_diag to get A*P */
1140   for (i = 0; i < am; i++) {
1141     aopnz  =  aopothi[i+1] -  aopothi[i];
1142     adponz = adpoi[i+1] - adpoi[i];
1143     adpdnz = adpdi[i+1] - adpdi[i];
1144 
1145     /* Correct indices from A_diag*P_diag */
1146     for (i1 = 0; i1 < adpdnz; i1++) {
1147       adpdJ[i1] += p_colstart;
1148     }
1149     /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1150     Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1151     ierr = MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz);CHKERRQ(ierr);
1152 
1153     aopJ += aopnz;
1154     adpoJ += adponz;
1155     adpdJ += adpdnz;
1156     apJ += apnz;
1157     api[i+1] = api[i] + apnz;
1158   }
1159 
1160   /* malloc apa to store dense row A[i,:]*P */
1161   ierr = PetscCalloc1(pN+2,&ptap->apa);CHKERRQ(ierr);
1162 
1163   /* create and assemble symbolic parallel matrix C */
1164   ierr = MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1165   ierr = MatSetBlockSizesFromMats(C,A,P);CHKERRQ(ierr);
1166   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1167   ierr = MatSetType(C,mtype);CHKERRQ(ierr);
1168   ierr = MatMPIAIJSetPreallocation(C,0,dnz,0,onz);CHKERRQ(ierr);
1169   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1170 
1171   ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);CHKERRQ(ierr);
1172   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1173   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1174   ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1175 
1176   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1177   C->ops->productnumeric = MatProductNumeric_AB;
1178 
1179   /* attach the supporting struct to C for reuse */
1180   C->product->data    = ptap;
1181   C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;
1182 
1183   /* set MatInfo */
1184   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
1185   if (afill < 1.0) afill = 1.0;
1186   C->info.mallocs           = nspacedouble;
1187   C->info.fill_ratio_given  = fill;
1188   C->info.fill_ratio_needed = afill;
1189 
1190 #if defined(PETSC_USE_INFO)
1191   if (api[am]) {
1192     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
1193     ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
1194   } else {
1195     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
1196   }
1197 #endif
1198 
1199   ierr = MatDestroy(&aopoth);CHKERRQ(ierr);
1200   ierr = MatDestroy(&adpd);CHKERRQ(ierr);
1201   ierr = PetscFree(j_temp);CHKERRQ(ierr);
1202   ierr = PetscFree(adpoj);CHKERRQ(ierr);
1203   ierr = PetscFree(adpoi);CHKERRQ(ierr);
1204   PetscFunctionReturn(0);
1205 }
1206 
1207 /*-------------------------------------------------------------------------*/
1208 /* This routine only works when scall=MAT_REUSE_MATRIX! */
1209 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
1210 {
1211   PetscErrorCode ierr;
1212   Mat_APMPI      *ptap;
1213   Mat            Pt;
1214 
1215   PetscFunctionBegin;
1216   MatCheckProduct(C,3);
1217   ptap = (Mat_APMPI*)C->product->data;
1218   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
1219   if (!ptap->Pt) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()");
1220 
1221   Pt   = ptap->Pt;
1222   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
1223   ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(Pt,A,C);CHKERRQ(ierr);
1224   PetscFunctionReturn(0);
1225 }
1226 
1227 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1228 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat C)
1229 {
1230   PetscErrorCode      ierr;
1231   Mat_APMPI           *ptap;
1232   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data;
1233   MPI_Comm            comm;
1234   PetscMPIInt         size,rank;
1235   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1236   PetscInt            pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n;
1237   PetscInt            *lnk,i,k,nsend,rstart;
1238   PetscBT             lnkbt;
1239   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,nrecv;
1240   PETSC_UNUSED PetscMPIInt icompleted=0;
1241   PetscInt            **buf_rj,**buf_ri,**buf_ri_k,row,ncols,*cols;
1242   PetscInt            len,proc,*dnz,*onz,*owners,nzi;
1243   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1244   MPI_Request         *swaits,*rwaits;
1245   MPI_Status          *sstatus,rstatus;
1246   PetscLayout         rowmap;
1247   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
1248   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
1249   PetscInt            *Jptr,*prmap=p->garray,con,j,Crmax;
1250   Mat_SeqAIJ          *a_loc,*c_loc,*c_oth;
1251   PetscTable          ta;
1252   MatType             mtype;
1253   const char          *prefix;
1254 
1255   PetscFunctionBegin;
1256   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1257   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
1258   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
1259 
1260   /* create symbolic parallel matrix C */
1261   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1262   ierr = MatSetType(C,mtype);CHKERRQ(ierr);
1263 
1264   C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1265 
1266   /* create struct Mat_APMPI and attached it to C later */
1267   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1268   ptap->reuse = MAT_INITIAL_MATRIX;
1269 
1270   /* (0) compute Rd = Pd^T, Ro = Po^T  */
1271   /* --------------------------------- */
1272   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1273   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1274 
1275   /* (1) compute symbolic A_loc */
1276   /* ---------------------------*/
1277   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);CHKERRQ(ierr);
1278 
1279   /* (2-1) compute symbolic C_oth = Ro*A_loc  */
1280   /* ------------------------------------ */
1281   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1282   ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr);
1283   ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr);
1284   ierr = MatCreate(PETSC_COMM_SELF,&ptap->C_oth);CHKERRQ(ierr);
1285   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,ptap->C_oth);CHKERRQ(ierr);
1286 
1287   /* (3) send coj of C_oth to other processors  */
1288   /* ------------------------------------------ */
1289   /* determine row ownership */
1290   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
1291   rowmap->n  = pn;
1292   rowmap->bs = 1;
1293   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
1294   owners = rowmap->range;
1295 
1296   /* determine the number of messages to send, their lengths */
1297   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
1298   ierr = PetscArrayzero(len_s,size);CHKERRQ(ierr);
1299   ierr = PetscArrayzero(len_si,size);CHKERRQ(ierr);
1300 
1301   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1302   coi   = c_oth->i; coj = c_oth->j;
1303   con   = ptap->C_oth->rmap->n;
1304   proc  = 0;
1305   for (i=0; i<con; i++) {
1306     while (prmap[i] >= owners[proc+1]) proc++;
1307     len_si[proc]++;               /* num of rows in Co(=Pt*A) to be sent to [proc] */
1308     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1309   }
1310 
1311   len          = 0; /* max length of buf_si[], see (4) */
1312   owners_co[0] = 0;
1313   nsend        = 0;
1314   for (proc=0; proc<size; proc++) {
1315     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1316     if (len_s[proc]) {
1317       nsend++;
1318       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1319       len         += len_si[proc];
1320     }
1321   }
1322 
1323   /* determine the number and length of messages to receive for coi and coj  */
1324   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
1325   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
1326 
1327   /* post the Irecv and Isend of coj */
1328   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1329   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1330   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
1331   for (proc=0, k=0; proc<size; proc++) {
1332     if (!len_s[proc]) continue;
1333     i    = owners_co[proc];
1334     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRMPI(ierr);
1335     k++;
1336   }
1337 
1338   /* (2-2) compute symbolic C_loc = Rd*A_loc */
1339   /* ---------------------------------------- */
1340   ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr);
1341   ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr);
1342   ierr = MatCreate(PETSC_COMM_SELF,&ptap->C_loc);CHKERRQ(ierr);
1343   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,ptap->C_loc);CHKERRQ(ierr);
1344   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
1345 
1346   /* receives coj are complete */
1347   for (i=0; i<nrecv; i++) {
1348     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr);
1349   }
1350   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1351   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRMPI(ierr);}
1352 
1353   /* add received column indices into ta to update Crmax */
1354   a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data;
1355 
1356   /* create and initialize a linked list */
1357   ierr = PetscTableCreate(an,aN,&ta);CHKERRQ(ierr); /* for compute Crmax */
1358   MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta);
1359 
1360   for (k=0; k<nrecv; k++) {/* k-th received message */
1361     Jptr = buf_rj[k];
1362     for (j=0; j<len_r[k]; j++) {
1363       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
1364     }
1365   }
1366   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
1367   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
1368 
1369   /* (4) send and recv coi */
1370   /*-----------------------*/
1371   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1372   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1373   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1374   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1375   for (proc=0,k=0; proc<size; proc++) {
1376     if (!len_s[proc]) continue;
1377     /* form outgoing message for i-structure:
1378          buf_si[0]:                 nrows to be sent
1379                [1:nrows]:           row index (global)
1380                [nrows+1:2*nrows+1]: i-structure index
1381     */
1382     /*-------------------------------------------*/
1383     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1384     buf_si_i    = buf_si + nrows+1;
1385     buf_si[0]   = nrows;
1386     buf_si_i[0] = 0;
1387     nrows       = 0;
1388     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1389       nzi = coi[i+1] - coi[i];
1390       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1391       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1392       nrows++;
1393     }
1394     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRMPI(ierr);
1395     k++;
1396     buf_si += len_si[proc];
1397   }
1398   for (i=0; i<nrecv; i++) {
1399     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr);
1400   }
1401   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1402   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRMPI(ierr);}
1403 
1404   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
1405   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1406   ierr = PetscFree(swaits);CHKERRQ(ierr);
1407   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1408 
1409   /* (5) compute the local portion of C      */
1410   /* ------------------------------------------ */
1411   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of C */
1412   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
1413   current_space = free_space;
1414 
1415   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
1416   for (k=0; k<nrecv; k++) {
1417     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1418     nrows       = *buf_ri_k[k];
1419     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1420     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1421   }
1422 
1423   ierr = MatPreallocateInitialize(comm,pn,an,dnz,onz);CHKERRQ(ierr);
1424   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1425   for (i=0; i<pn; i++) { /* for each local row of C */
1426     /* add C_loc into C */
1427     nzi  = c_loc->i[i+1] - c_loc->i[i];
1428     Jptr = c_loc->j + c_loc->i[i];
1429     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1430 
1431     /* add received col data into lnk */
1432     for (k=0; k<nrecv; k++) { /* k-th received message */
1433       if (i == *nextrow[k]) { /* i-th row */
1434         nzi  = *(nextci[k]+1) - *nextci[k];
1435         Jptr = buf_rj[k] + *nextci[k];
1436         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1437         nextrow[k]++; nextci[k]++;
1438       }
1439     }
1440 
1441     /* add missing diagonal entry */
1442     if (C->force_diagonals) {
1443       k = i + owners[rank]; /* column index */
1444       ierr = PetscLLCondensedAddSorted(1,&k,lnk,lnkbt);CHKERRQ(ierr);
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(C,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1459   if (P->cmap->bs > 0) {ierr = PetscLayoutSetBlockSize(C->rmap,P->cmap->bs);CHKERRQ(ierr);}
1460   if (A->cmap->bs > 0) {ierr = PetscLayoutSetBlockSize(C->cmap,A->cmap->bs);CHKERRQ(ierr);}
1461   ierr = MatMPIAIJSetPreallocation(C,0,dnz,0,onz);CHKERRQ(ierr);
1462   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1463 
1464   /* add C_loc and C_oth to C */
1465   ierr = MatGetOwnershipRange(C,&rstart,NULL);CHKERRQ(ierr);
1466   for (i=0; i<pn; i++) {
1467     ncols = c_loc->i[i+1] - c_loc->i[i];
1468     cols  = c_loc->j + c_loc->i[i];
1469     row   = rstart + i;
1470     ierr = MatSetValues(C,1,(const PetscInt*)&row,ncols,(const PetscInt*)cols,NULL,INSERT_VALUES);CHKERRQ(ierr);
1471 
1472     if (C->force_diagonals) {
1473       ierr = MatSetValues(C,1,(const PetscInt*)&row,1,(const PetscInt*)&row,NULL,INSERT_VALUES);CHKERRQ(ierr);
1474     }
1475   }
1476   for (i=0; i<con; i++) {
1477     ncols = c_oth->i[i+1] - c_oth->i[i];
1478     cols  = c_oth->j + c_oth->i[i];
1479     row   = prmap[i];
1480     ierr = MatSetValues(C,1,(const PetscInt*)&row,ncols,(const PetscInt*)cols,NULL,INSERT_VALUES);CHKERRQ(ierr);
1481   }
1482   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1483   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1484   ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1485 
1486   /* members in merge */
1487   ierr = PetscFree(id_r);CHKERRQ(ierr);
1488   ierr = PetscFree(len_r);CHKERRQ(ierr);
1489   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
1490   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
1491   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
1492   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
1493   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
1494 
1495   /* attach the supporting struct to C for reuse */
1496   C->product->data    = ptap;
1497   C->product->destroy = MatDestroy_MPIAIJ_PtAP;
1498   PetscFunctionReturn(0);
1499 }
1500 
1501 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
1502 {
1503   PetscErrorCode    ierr;
1504   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data;
1505   Mat_SeqAIJ        *c_seq;
1506   Mat_APMPI         *ptap;
1507   Mat               A_loc,C_loc,C_oth;
1508   PetscInt          i,rstart,rend,cm,ncols,row;
1509   const PetscInt    *cols;
1510   const PetscScalar *vals;
1511 
1512   PetscFunctionBegin;
1513   MatCheckProduct(C,3);
1514   ptap = (Mat_APMPI*)C->product->data;
1515   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
1516   if (!ptap->A_loc) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()");
1517   ierr = MatZeroEntries(C);CHKERRQ(ierr);
1518 
1519   if (ptap->reuse == MAT_REUSE_MATRIX) {
1520     /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1521     /* 1) get R = Pd^T, Ro = Po^T */
1522     /*----------------------------*/
1523     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1524     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1525 
1526     /* 2) compute numeric A_loc */
1527     /*--------------------------*/
1528     ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);CHKERRQ(ierr);
1529   }
1530 
1531   /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
1532   A_loc = ptap->A_loc;
1533   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);CHKERRQ(ierr);
1534   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);CHKERRQ(ierr);
1535   C_loc = ptap->C_loc;
1536   C_oth = ptap->C_oth;
1537 
1538   /* add C_loc and C_oth to C */
1539   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
1540 
1541   /* C_loc -> C */
1542   cm    = C_loc->rmap->N;
1543   c_seq = (Mat_SeqAIJ*)C_loc->data;
1544   cols = c_seq->j;
1545   vals = c_seq->a;
1546   for (i=0; i<cm; i++) {
1547     ncols = c_seq->i[i+1] - c_seq->i[i];
1548     row = rstart + i;
1549     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1550     cols += ncols; vals += ncols;
1551   }
1552 
1553   /* Co -> C, off-processor part */
1554   cm    = C_oth->rmap->N;
1555   c_seq = (Mat_SeqAIJ*)C_oth->data;
1556   cols  = c_seq->j;
1557   vals  = c_seq->a;
1558   for (i=0; i<cm; i++) {
1559     ncols = c_seq->i[i+1] - c_seq->i[i];
1560     row = p->garray[i];
1561     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1562     cols += ncols; vals += ncols;
1563   }
1564   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1565   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1566   ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1567 
1568   ptap->reuse = MAT_REUSE_MATRIX;
1569   PetscFunctionReturn(0);
1570 }
1571 
1572 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1573 {
1574   PetscErrorCode      ierr;
1575   Mat_Merge_SeqsToMPI *merge;
1576   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data;
1577   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1578   Mat_APMPI           *ptap;
1579   PetscInt            *adj;
1580   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1581   MatScalar           *ada,*ca,valtmp;
1582   PetscInt            am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1583   MPI_Comm            comm;
1584   PetscMPIInt         size,rank,taga,*len_s;
1585   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1586   PetscInt            **buf_ri,**buf_rj;
1587   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1588   MPI_Request         *s_waits,*r_waits;
1589   MPI_Status          *status;
1590   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1591   PetscInt            *ai,*aj,*coi,*coj,*poJ,*pdJ;
1592   Mat                 A_loc;
1593   Mat_SeqAIJ          *a_loc;
1594 
1595   PetscFunctionBegin;
1596   MatCheckProduct(C,3);
1597   ptap = (Mat_APMPI*)C->product->data;
1598   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
1599   if (!ptap->A_loc) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()");
1600   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1601   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
1602   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
1603 
1604   merge = ptap->merge;
1605 
1606   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1607   /*------------------------------------------*/
1608   /* get data from symbolic products */
1609   coi    = merge->coi; coj = merge->coj;
1610   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1611   bi     = merge->bi; bj = merge->bj;
1612   owners = merge->rowmap->range;
1613   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1614 
1615   /* get A_loc by taking all local rows of A */
1616   A_loc = ptap->A_loc;
1617   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1618   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1619   ai    = a_loc->i;
1620   aj    = a_loc->j;
1621 
1622   for (i=0; i<am; i++) {
1623     anz = ai[i+1] - ai[i];
1624     adj = aj + ai[i];
1625     ada = a_loc->a + ai[i];
1626 
1627     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1628     /*-------------------------------------------------------------*/
1629     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1630     pnz = po->i[i+1] - po->i[i];
1631     poJ = po->j + po->i[i];
1632     pA  = po->a + po->i[i];
1633     for (j=0; j<pnz; j++) {
1634       row = poJ[j];
1635       cj  = coj + coi[row];
1636       ca  = coa + coi[row];
1637       /* perform sparse axpy */
1638       nexta  = 0;
1639       valtmp = pA[j];
1640       for (k=0; nexta<anz; k++) {
1641         if (cj[k] == adj[nexta]) {
1642           ca[k] += valtmp*ada[nexta];
1643           nexta++;
1644         }
1645       }
1646       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1647     }
1648 
1649     /* put the value into Cd (diagonal part) */
1650     pnz = pd->i[i+1] - pd->i[i];
1651     pdJ = pd->j + pd->i[i];
1652     pA  = pd->a + pd->i[i];
1653     for (j=0; j<pnz; j++) {
1654       row = pdJ[j];
1655       cj  = bj + bi[row];
1656       ca  = ba + bi[row];
1657       /* perform sparse axpy */
1658       nexta  = 0;
1659       valtmp = pA[j];
1660       for (k=0; nexta<anz; k++) {
1661         if (cj[k] == adj[nexta]) {
1662           ca[k] += valtmp*ada[nexta];
1663           nexta++;
1664         }
1665       }
1666       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1667     }
1668   }
1669 
1670   /* 3) send and recv matrix values coa */
1671   /*------------------------------------*/
1672   buf_ri = merge->buf_ri;
1673   buf_rj = merge->buf_rj;
1674   len_s  = merge->len_s;
1675   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1676   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1677 
1678   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1679   for (proc=0,k=0; proc<size; proc++) {
1680     if (!len_s[proc]) continue;
1681     i    = merge->owners_co[proc];
1682     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRMPI(ierr);
1683     k++;
1684   }
1685   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRMPI(ierr);}
1686   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRMPI(ierr);}
1687 
1688   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1689   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1690   ierr = PetscFree(coa);CHKERRQ(ierr);
1691 
1692   /* 4) insert local Cseq and received values into Cmpi */
1693   /*----------------------------------------------------*/
1694   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1695   for (k=0; k<merge->nrecv; k++) {
1696     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1697     nrows       = *(buf_ri_k[k]);
1698     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1699     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1700   }
1701 
1702   for (i=0; i<cm; i++) {
1703     row  = owners[rank] + i; /* global row index of C_seq */
1704     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1705     ba_i = ba + bi[i];
1706     bnz  = bi[i+1] - bi[i];
1707     /* add received vals into ba */
1708     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1709       /* i-th row */
1710       if (i == *nextrow[k]) {
1711         cnz    = *(nextci[k]+1) - *nextci[k];
1712         cj     = buf_rj[k] + *(nextci[k]);
1713         ca     = abuf_r[k] + *(nextci[k]);
1714         nextcj = 0;
1715         for (j=0; nextcj<cnz; j++) {
1716           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1717             ba_i[j] += ca[nextcj++];
1718           }
1719         }
1720         nextrow[k]++; nextci[k]++;
1721         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1722       }
1723     }
1724     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1725   }
1726   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1727   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1728 
1729   ierr = PetscFree(ba);CHKERRQ(ierr);
1730   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1731   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1732   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1733   PetscFunctionReturn(0);
1734 }
1735 
1736 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat C)
1737 {
1738   PetscErrorCode      ierr;
1739   Mat                 A_loc,POt,PDt;
1740   Mat_APMPI           *ptap;
1741   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1742   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data;
1743   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1744   PetscInt            nnz;
1745   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1746   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1747   MPI_Comm            comm;
1748   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1749   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1750   PetscInt            len,proc,*dnz,*onz,*owners;
1751   PetscInt            nzi,*bi,*bj;
1752   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1753   MPI_Request         *swaits,*rwaits;
1754   MPI_Status          *sstatus,rstatus;
1755   Mat_Merge_SeqsToMPI *merge;
1756   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1757   PetscReal           afill  =1.0,afill_tmp;
1758   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1759   Mat_SeqAIJ          *a_loc,*pdt,*pot;
1760   PetscTable          ta;
1761   MatType             mtype;
1762 
1763   PetscFunctionBegin;
1764   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1765   /* check if matrix local sizes are compatible */
1766   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);
1767 
1768   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
1769   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
1770 
1771   /* create struct Mat_APMPI and attached it to C later */
1772   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1773 
1774   /* get A_loc by taking all local rows of A */
1775   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1776 
1777   ptap->A_loc = A_loc;
1778   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1779   ai          = a_loc->i;
1780   aj          = a_loc->j;
1781 
1782   /* determine symbolic Co=(p->B)^T*A - send to others */
1783   /*----------------------------------------------------*/
1784   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1785   pdt  = (Mat_SeqAIJ*)PDt->data;
1786   pdti = pdt->i; pdtj = pdt->j;
1787 
1788   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1789   pot  = (Mat_SeqAIJ*)POt->data;
1790   poti = pot->i; potj = pot->j;
1791 
1792   /* then, compute symbolic Co = (p->B)^T*A */
1793   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1794                          >= (num of nonzero rows of C_seq) - pn */
1795   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1796   coi[0] = 0;
1797 
1798   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1799   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1800   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1801   current_space = free_space;
1802 
1803   /* create and initialize a linked list */
1804   ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr);
1805   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
1806   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
1807 
1808   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1809 
1810   for (i=0; i<pon; i++) {
1811     pnz = poti[i+1] - poti[i];
1812     ptJ = potj + poti[i];
1813     for (j=0; j<pnz; j++) {
1814       row  = ptJ[j]; /* row of A_loc == col of Pot */
1815       anz  = ai[row+1] - ai[row];
1816       Jptr = aj + ai[row];
1817       /* add non-zero cols of AP into the sorted linked list lnk */
1818       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1819     }
1820     nnz = lnk[0];
1821 
1822     /* If free space is not available, double the total space in the list */
1823     if (current_space->local_remaining<nnz) {
1824       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1825       nspacedouble++;
1826     }
1827 
1828     /* Copy data into free space, and zero out denserows */
1829     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1830 
1831     current_space->array           += nnz;
1832     current_space->local_used      += nnz;
1833     current_space->local_remaining -= nnz;
1834 
1835     coi[i+1] = coi[i] + nnz;
1836   }
1837 
1838   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1839   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1840   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
1841 
1842   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1843   if (afill_tmp > afill) afill = afill_tmp;
1844 
1845   /* send j-array (coj) of Co to other processors */
1846   /*----------------------------------------------*/
1847   /* determine row ownership */
1848   ierr = PetscNew(&merge);CHKERRQ(ierr);
1849   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1850 
1851   merge->rowmap->n  = pn;
1852   merge->rowmap->bs = 1;
1853 
1854   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1855   owners = merge->rowmap->range;
1856 
1857   /* determine the number of messages to send, their lengths */
1858   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1859   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
1860 
1861   len_s        = merge->len_s;
1862   merge->nsend = 0;
1863 
1864   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1865 
1866   proc = 0;
1867   for (i=0; i<pon; i++) {
1868     while (prmap[i] >= owners[proc+1]) proc++;
1869     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1870     len_s[proc] += coi[i+1] - coi[i];
1871   }
1872 
1873   len          = 0; /* max length of buf_si[] */
1874   owners_co[0] = 0;
1875   for (proc=0; proc<size; proc++) {
1876     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1877     if (len_si[proc]) {
1878       merge->nsend++;
1879       len_si[proc] = 2*(len_si[proc] + 1);
1880       len         += len_si[proc];
1881     }
1882   }
1883 
1884   /* determine the number and length of messages to receive for coi and coj  */
1885   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1886   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1887 
1888   /* post the Irecv and Isend of coj */
1889   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1890   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1891   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1892   for (proc=0, k=0; proc<size; proc++) {
1893     if (!len_s[proc]) continue;
1894     i    = owners_co[proc];
1895     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRMPI(ierr);
1896     k++;
1897   }
1898 
1899   /* receives and sends of coj are complete */
1900   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1901   for (i=0; i<merge->nrecv; i++) {
1902     PETSC_UNUSED PetscMPIInt icompleted;
1903     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr);
1904   }
1905   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1906   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRMPI(ierr);}
1907 
1908   /* add received column indices into table to update Armax */
1909   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/tutorials/ex56.c! */
1910   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
1911     Jptr = buf_rj[k];
1912     for (j=0; j<merge->len_r[k]; j++) {
1913       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
1914     }
1915   }
1916   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
1917   /* 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); */
1918 
1919   /* send and recv coi */
1920   /*-------------------*/
1921   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1922   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1923   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1924   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1925   for (proc=0,k=0; proc<size; proc++) {
1926     if (!len_s[proc]) continue;
1927     /* form outgoing message for i-structure:
1928          buf_si[0]:                 nrows to be sent
1929                [1:nrows]:           row index (global)
1930                [nrows+1:2*nrows+1]: i-structure index
1931     */
1932     /*-------------------------------------------*/
1933     nrows       = len_si[proc]/2 - 1;
1934     buf_si_i    = buf_si + nrows+1;
1935     buf_si[0]   = nrows;
1936     buf_si_i[0] = 0;
1937     nrows       = 0;
1938     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1939       nzi               = coi[i+1] - coi[i];
1940       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1941       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1942       nrows++;
1943     }
1944     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRMPI(ierr);
1945     k++;
1946     buf_si += len_si[proc];
1947   }
1948   i = merge->nrecv;
1949   while (i--) {
1950     PETSC_UNUSED PetscMPIInt icompleted;
1951     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr);
1952   }
1953   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1954   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRMPI(ierr);}
1955   ierr = PetscFree(len_si);CHKERRQ(ierr);
1956   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1957   ierr = PetscFree(swaits);CHKERRQ(ierr);
1958   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1959   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1960 
1961   /* compute the local portion of C (mpi mat) */
1962   /*------------------------------------------*/
1963   /* allocate bi array and free space for accumulating nonzero column info */
1964   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1965   bi[0] = 0;
1966 
1967   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1968   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1969   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1970   current_space = free_space;
1971 
1972   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1973   for (k=0; k<merge->nrecv; k++) {
1974     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1975     nrows       = *buf_ri_k[k];
1976     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1977     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th received i-structure  */
1978   }
1979 
1980   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1981   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1982   rmax = 0;
1983   for (i=0; i<pn; i++) {
1984     /* add pdt[i,:]*AP into lnk */
1985     pnz = pdti[i+1] - pdti[i];
1986     ptJ = pdtj + pdti[i];
1987     for (j=0; j<pnz; j++) {
1988       row  = ptJ[j];  /* row of AP == col of Pt */
1989       anz  = ai[row+1] - ai[row];
1990       Jptr = aj + ai[row];
1991       /* add non-zero cols of AP into the sorted linked list lnk */
1992       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1993     }
1994 
1995     /* add received col data into lnk */
1996     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1997       if (i == *nextrow[k]) { /* i-th row */
1998         nzi  = *(nextci[k]+1) - *nextci[k];
1999         Jptr = buf_rj[k] + *nextci[k];
2000         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
2001         nextrow[k]++; nextci[k]++;
2002       }
2003     }
2004 
2005     /* add missing diagonal entry */
2006     if (C->force_diagonals) {
2007       k = i + owners[rank]; /* column index */
2008       ierr = PetscLLCondensedAddSorted_Scalable(1,&k,lnk);CHKERRQ(ierr);
2009     }
2010 
2011     nnz = lnk[0];
2012 
2013     /* if free space is not available, make more free space */
2014     if (current_space->local_remaining<nnz) {
2015       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
2016       nspacedouble++;
2017     }
2018     /* copy data into free space, then initialize lnk */
2019     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
2020     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2021 
2022     current_space->array           += nnz;
2023     current_space->local_used      += nnz;
2024     current_space->local_remaining -= nnz;
2025 
2026     bi[i+1] = bi[i] + nnz;
2027     if (nnz > rmax) rmax = nnz;
2028   }
2029   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
2030 
2031   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
2032   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
2033   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
2034   if (afill_tmp > afill) afill = afill_tmp;
2035   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
2036   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
2037 
2038   ierr = MatDestroy(&POt);CHKERRQ(ierr);
2039   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
2040 
2041   /* create symbolic parallel matrix C - why cannot be assembled in Numeric part   */
2042   /*-------------------------------------------------------------------------------*/
2043   ierr = MatSetSizes(C,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
2044   ierr = MatSetBlockSizes(C,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
2045   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
2046   ierr = MatSetType(C,mtype);CHKERRQ(ierr);
2047   ierr = MatMPIAIJSetPreallocation(C,0,dnz,0,onz);CHKERRQ(ierr);
2048   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2049   ierr = MatSetBlockSize(C,1);CHKERRQ(ierr);
2050   ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
2051   for (i=0; i<pn; i++) {
2052     row  = i + rstart;
2053     nnz  = bi[i+1] - bi[i];
2054     Jptr = bj + bi[i];
2055     ierr = MatSetValues(C,1,&row,nnz,Jptr,NULL,INSERT_VALUES);CHKERRQ(ierr);
2056   }
2057   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2058   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2059   ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
2060   merge->bi        = bi;
2061   merge->bj        = bj;
2062   merge->coi       = coi;
2063   merge->coj       = coj;
2064   merge->buf_ri    = buf_ri;
2065   merge->buf_rj    = buf_rj;
2066   merge->owners_co = owners_co;
2067 
2068   /* attach the supporting struct to C for reuse */
2069   C->product->data    = ptap;
2070   C->product->destroy = MatDestroy_MPIAIJ_PtAP;
2071   ptap->merge         = merge;
2072 
2073   C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
2074 
2075 #if defined(PETSC_USE_INFO)
2076   if (bi[pn] != 0) {
2077     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
2078     ierr = PetscInfo1(C,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
2079   } else {
2080     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
2081   }
2082 #endif
2083   PetscFunctionReturn(0);
2084 }
2085 
2086 /* ---------------------------------------------------------------- */
2087 static PetscErrorCode MatProductSymbolic_AtB_MPIAIJ_MPIAIJ(Mat C)
2088 {
2089   PetscErrorCode ierr;
2090   Mat_Product    *product = C->product;
2091   Mat            A=product->A,B=product->B;
2092   PetscReal      fill=product->fill;
2093   PetscBool      flg;
2094 
2095   PetscFunctionBegin;
2096   /* scalable */
2097   ierr = PetscStrcmp(product->alg,"scalable",&flg);CHKERRQ(ierr);
2098   if (flg) {
2099     ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
2100     goto next;
2101   }
2102 
2103   /* nonscalable */
2104   ierr = PetscStrcmp(product->alg,"nonscalable",&flg);CHKERRQ(ierr);
2105   if (flg) {
2106     ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
2107     goto next;
2108   }
2109 
2110   /* matmatmult */
2111   ierr = PetscStrcmp(product->alg,"at*b",&flg);CHKERRQ(ierr);
2112   if (flg) {
2113     Mat       At;
2114     Mat_APMPI *ptap;
2115 
2116     ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr);
2117     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(At,B,fill,C);CHKERRQ(ierr);
2118     ptap = (Mat_APMPI*)C->product->data;
2119     if (ptap) {
2120       ptap->Pt = At;
2121       C->product->destroy = MatDestroy_MPIAIJ_PtAP;
2122     }
2123     C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
2124     goto next;
2125   }
2126 
2127   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProduct type is not supported");
2128 
2129 next:
2130   C->ops->productnumeric = MatProductNumeric_AtB;
2131   PetscFunctionReturn(0);
2132 }
2133 
2134 /* ---------------------------------------------------------------- */
2135 /* Set options for MatMatMultxxx_MPIAIJ_MPIAIJ */
2136 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AB(Mat C)
2137 {
2138   PetscErrorCode ierr;
2139   Mat_Product    *product = C->product;
2140   Mat            A=product->A,B=product->B;
2141 #if defined(PETSC_HAVE_HYPRE)
2142   const char     *algTypes[4] = {"scalable","nonscalable","seqmpi","hypre"};
2143   PetscInt       nalg = 4;
2144 #else
2145   const char     *algTypes[3] = {"scalable","nonscalable","seqmpi"};
2146   PetscInt       nalg = 3;
2147 #endif
2148   PetscInt       alg = 1; /* set nonscalable algorithm as default */
2149   PetscBool      flg;
2150   MPI_Comm       comm;
2151 
2152   PetscFunctionBegin;
2153   /* Check matrix local sizes */
2154   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
2155   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(PETSC_COMM_SELF,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);
2156 
2157   /* Set "nonscalable" as default algorithm */
2158   ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr);
2159   if (flg) {
2160     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2161 
2162     /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2163     if (B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
2164       MatInfo     Ainfo,Binfo;
2165       PetscInt    nz_local;
2166       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
2167 
2168       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
2169       ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr);
2170       nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);
2171 
2172       if (B->cmap->N > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE;
2173       ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
2174 
2175       if (alg_scalable) {
2176         alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2177         ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2178         ierr = PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,product->fill*nz_local);CHKERRQ(ierr);
2179       }
2180     }
2181   }
2182 
2183   /* Get runtime option */
2184   if (product->api_user) {
2185     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMult","Mat");CHKERRQ(ierr);
2186     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2187     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2188   } else {
2189     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AB","Mat");CHKERRQ(ierr);
2190     ierr = PetscOptionsEList("-matproduct_ab_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2191     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2192   }
2193   if (flg) {
2194     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2195   }
2196 
2197   C->ops->productsymbolic = MatProductSymbolic_AB_MPIAIJ_MPIAIJ;
2198   PetscFunctionReturn(0);
2199 }
2200 
2201 /* Set options for MatTransposeMatMultXXX_MPIAIJ_MPIAIJ */
2202 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AtB(Mat C)
2203 {
2204   PetscErrorCode ierr;
2205   Mat_Product    *product = C->product;
2206   Mat            A=product->A,B=product->B;
2207   const char     *algTypes[3] = {"scalable","nonscalable","at*b"};
2208   PetscInt       nalg = 3;
2209   PetscInt       alg = 1; /* set default algorithm  */
2210   PetscBool      flg;
2211   MPI_Comm       comm;
2212 
2213   PetscFunctionBegin;
2214   /* Check matrix local sizes */
2215   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
2216   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != B (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);
2217 
2218   /* Set default algorithm */
2219   ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr);
2220   if (flg) {
2221     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2222   }
2223 
2224   /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2225   if (alg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
2226     MatInfo     Ainfo,Binfo;
2227     PetscInt    nz_local;
2228     PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
2229 
2230     ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
2231     ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr);
2232     nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);
2233 
2234     if (B->cmap->N > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE;
2235     ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
2236 
2237     if (alg_scalable) {
2238       alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2239       ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2240       ierr = PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,product->fill*nz_local);CHKERRQ(ierr);
2241     }
2242   }
2243 
2244   /* Get runtime option */
2245   if (product->api_user) {
2246     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatTransposeMatMult","Mat");CHKERRQ(ierr);
2247     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2248     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2249   } else {
2250     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AtB","Mat");CHKERRQ(ierr);
2251     ierr = PetscOptionsEList("-matproduct_atb_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2252     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2253   }
2254   if (flg) {
2255     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2256   }
2257 
2258   C->ops->productsymbolic = MatProductSymbolic_AtB_MPIAIJ_MPIAIJ;
2259   PetscFunctionReturn(0);
2260 }
2261 
2262 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_PtAP(Mat C)
2263 {
2264   PetscErrorCode ierr;
2265   Mat_Product    *product = C->product;
2266   Mat            A=product->A,P=product->B;
2267   MPI_Comm       comm;
2268   PetscBool      flg;
2269   PetscInt       alg=1; /* set default algorithm */
2270 #if !defined(PETSC_HAVE_HYPRE)
2271   const char     *algTypes[4] = {"scalable","nonscalable","allatonce","allatonce_merged"};
2272   PetscInt       nalg=4;
2273 #else
2274   const char     *algTypes[5] = {"scalable","nonscalable","allatonce","allatonce_merged","hypre"};
2275   PetscInt       nalg=5;
2276 #endif
2277   PetscInt       pN=P->cmap->N;
2278 
2279   PetscFunctionBegin;
2280   /* Check matrix local sizes */
2281   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
2282   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
2283   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
2284 
2285   /* Set "nonscalable" as default algorithm */
2286   ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr);
2287   if (flg) {
2288     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2289 
2290     /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2291     if (pN > 100000) {
2292       MatInfo     Ainfo,Pinfo;
2293       PetscInt    nz_local;
2294       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
2295 
2296       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
2297       ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
2298       nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
2299 
2300       if (pN > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE;
2301       ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
2302 
2303       if (alg_scalable) {
2304         alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2305         ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2306       }
2307     }
2308   }
2309 
2310   /* Get runtime option */
2311   if (product->api_user) {
2312     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
2313     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2314     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2315   } else {
2316     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_PtAP","Mat");CHKERRQ(ierr);
2317     ierr = PetscOptionsEList("-matproduct_ptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2318     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2319   }
2320   if (flg) {
2321     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2322   }
2323 
2324   C->ops->productsymbolic = MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ;
2325   PetscFunctionReturn(0);
2326 }
2327 
2328 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_RARt(Mat C)
2329 {
2330   Mat_Product *product = C->product;
2331   Mat         A = product->A,R=product->B;
2332 
2333   PetscFunctionBegin;
2334   /* Check matrix local sizes */
2335   if (A->cmap->n != R->cmap->n || A->rmap->n != R->cmap->n) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A local (%D, %D), R local (%D,%D)",A->rmap->n,A->rmap->n,R->rmap->n,R->cmap->n);
2336 
2337   C->ops->productsymbolic = MatProductSymbolic_RARt_MPIAIJ_MPIAIJ;
2338   PetscFunctionReturn(0);
2339 }
2340 
2341 /*
2342  Set options for ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm
2343 */
2344 static PetscErrorCode MatProductSetFromOptions_MPIAIJ_ABC(Mat C)
2345 {
2346   PetscErrorCode ierr;
2347   Mat_Product    *product = C->product;
2348   PetscBool      flg = PETSC_FALSE;
2349   PetscInt       alg = 1; /* default algorithm */
2350   const char     *algTypes[3] = {"scalable","nonscalable","seqmpi"};
2351   PetscInt       nalg = 3;
2352 
2353   PetscFunctionBegin;
2354   /* Set default algorithm */
2355   ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr);
2356   if (flg) {
2357     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2358   }
2359 
2360   /* Get runtime option */
2361   if (product->api_user) {
2362     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMatMult","Mat");CHKERRQ(ierr);
2363     ierr = PetscOptionsEList("-matmatmatmult_via","Algorithmic approach","MatMatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2364     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2365   } else {
2366     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABC","Mat");CHKERRQ(ierr);
2367     ierr = PetscOptionsEList("-matproduct_abc_via","Algorithmic approach","MatProduct_ABC",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2368     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2369   }
2370   if (flg) {
2371     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2372   }
2373 
2374   C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ;
2375   C->ops->productsymbolic    = MatProductSymbolic_ABC;
2376   PetscFunctionReturn(0);
2377 }
2378 
2379 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat C)
2380 {
2381   PetscErrorCode ierr;
2382   Mat_Product    *product = C->product;
2383 
2384   PetscFunctionBegin;
2385   switch (product->type) {
2386   case MATPRODUCT_AB:
2387     ierr = MatProductSetFromOptions_MPIAIJ_AB(C);CHKERRQ(ierr);
2388     break;
2389   case MATPRODUCT_AtB:
2390     ierr = MatProductSetFromOptions_MPIAIJ_AtB(C);CHKERRQ(ierr);
2391     break;
2392   case MATPRODUCT_PtAP:
2393     ierr = MatProductSetFromOptions_MPIAIJ_PtAP(C);CHKERRQ(ierr);
2394     break;
2395   case MATPRODUCT_RARt:
2396     ierr = MatProductSetFromOptions_MPIAIJ_RARt(C);CHKERRQ(ierr);
2397     break;
2398   case MATPRODUCT_ABC:
2399     ierr = MatProductSetFromOptions_MPIAIJ_ABC(C);CHKERRQ(ierr);
2400     break;
2401   default:
2402     break;
2403   }
2404   PetscFunctionReturn(0);
2405 }
2406