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