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