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