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