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