xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision a7fac7c2b47dbcd8ece0cc2dfdfe0e63be1bb7b5)
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 
13 #undef __FUNCT__
14 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
15 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
16 {
17   PetscErrorCode ierr;
18   const char     *algTypes[2] = {"scalable","nonscalable"};
19   PetscInt       alg=0; /* set default algorithm */
20 
21   PetscFunctionBegin;
22   if (scall == MAT_INITIAL_MATRIX) {
23     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
24     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,2,algTypes[0],&alg,NULL);CHKERRQ(ierr);
25     ierr = PetscOptionsEnd();CHKERRQ(ierr);
26 
27     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
28     switch (alg) {
29     case 1:
30       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
31       break;
32     default:
33       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
34       break;
35     }
36     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
37   }
38   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
39   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
40   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
41   PetscFunctionReturn(0);
42 }
43 
44 #undef __FUNCT__
45 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
46 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
47 {
48   PetscErrorCode ierr;
49   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
50   Mat_PtAPMPI    *ptap = a->ptap;
51 
52   PetscFunctionBegin;
53   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
54   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
55   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
56   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
57   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
58   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
59   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
60   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
61   ierr = ptap->destroy(A);CHKERRQ(ierr);
62   ierr = PetscFree(ptap);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
68 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
69 {
70   PetscErrorCode ierr;
71   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
72   Mat_PtAPMPI    *ptap = a->ptap;
73 
74   PetscFunctionBegin;
75   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
76 
77   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
78   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
79   PetscFunctionReturn(0);
80 }
81 
82 #undef __FUNCT__
83 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
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   PetscInt       *adi=ad->i,*adj,*aoi=ao->i,*aoj;
91   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
92   Mat_SeqAIJ     *p_loc,*p_oth;
93   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
94   PetscScalar    *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
95   PetscInt       cm   =C->rmap->n,anz,pnz;
96   Mat_PtAPMPI    *ptap=c->ptap;
97   PetscInt       *api,*apj,*apJ,i,j,k,row;
98   PetscInt       cstart=C->cmap->rstart;
99   PetscInt       cdnz,conz,k0,k1;
100   MPI_Comm       comm;
101   PetscMPIInt    size;
102 
103   PetscFunctionBegin;
104   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
105   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
106 
107   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
108   /*-----------------------------------------------------*/
109   /* update numerical values of P_oth and P_loc */
110   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
111   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
112 
113   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
114   /*----------------------------------------------------------*/
115   /* get data from symbolic products */
116   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
117   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
118   if (size >1) {
119     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
120     pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
121   } else {
122     pi_oth=NULL; pj_oth=NULL; pa_oth=NULL;
123   }
124 
125   /* get apa for storing dense row A[i,:]*P */
126   apa = ptap->apa;
127 
128   api = ptap->api;
129   apj = ptap->apj;
130   for (i=0; i<cm; i++) {
131     /* diagonal portion of A */
132     anz = adi[i+1] - adi[i];
133     adj = ad->j + adi[i];
134     ada = ad->a + adi[i];
135     for (j=0; j<anz; j++) {
136       row = adj[j];
137       pnz = pi_loc[row+1] - pi_loc[row];
138       pj  = pj_loc + pi_loc[row];
139       pa  = pa_loc + pi_loc[row];
140 
141       /* perform dense axpy */
142       valtmp = ada[j];
143       for (k=0; k<pnz; k++) {
144         apa[pj[k]] += valtmp*pa[k];
145       }
146       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
147     }
148 
149     /* off-diagonal portion of A */
150     anz = aoi[i+1] - aoi[i];
151     aoj = ao->j + aoi[i];
152     aoa = ao->a + aoi[i];
153     for (j=0; j<anz; j++) {
154       row = aoj[j];
155       pnz = pi_oth[row+1] - pi_oth[row];
156       pj  = pj_oth + pi_oth[row];
157       pa  = pa_oth + pi_oth[row];
158 
159       /* perform dense axpy */
160       valtmp = aoa[j];
161       for (k=0; k<pnz; k++) {
162         apa[pj[k]] += valtmp*pa[k];
163       }
164       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
165     }
166 
167     /* set values in C */
168     apJ  = apj + api[i];
169     cdnz = cd->i[i+1] - cd->i[i];
170     conz = co->i[i+1] - co->i[i];
171 
172     /* 1st off-diagoanl part of C */
173     ca = coa + co->i[i];
174     k  = 0;
175     for (k0=0; k0<conz; k0++) {
176       if (apJ[k] >= cstart) break;
177       ca[k0]      = apa[apJ[k]];
178       apa[apJ[k]] = 0.0;
179       k++;
180     }
181 
182     /* diagonal part of C */
183     ca = cda + cd->i[i];
184     for (k1=0; k1<cdnz; k1++) {
185       ca[k1]      = apa[apJ[k]];
186       apa[apJ[k]] = 0.0;
187       k++;
188     }
189 
190     /* 2nd off-diagoanl part of C */
191     ca = coa + co->i[i];
192     for (; k0<conz; k0++) {
193       ca[k0]      = apa[apJ[k]];
194       apa[apJ[k]] = 0.0;
195       k++;
196     }
197   }
198   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
199   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
205 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
206 {
207   PetscErrorCode     ierr;
208   MPI_Comm           comm;
209   PetscMPIInt        size;
210   Mat                Cmpi;
211   Mat_PtAPMPI        *ptap;
212   PetscFreeSpaceList free_space=NULL,current_space=NULL;
213   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
214   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
215   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
216   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
217   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
218   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
219   PetscBT            lnkbt;
220   PetscScalar        *apa;
221   PetscReal          afill;
222   PetscInt           nlnk_max,armax,prmax;
223 
224   PetscFunctionBegin;
225   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
226   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
227 
228   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
229     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
230   }
231 
232   /* create struct Mat_PtAPMPI and attached it to C later */
233   ierr = PetscNew(&ptap);CHKERRQ(ierr);
234 
235   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
236   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
237 
238   /* get P_loc by taking all local rows of P */
239   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
240 
241   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
242   pi_loc = p_loc->i; pj_loc = p_loc->j;
243   if (size > 1) {
244     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
245     pi_oth = p_oth->i; pj_oth = p_oth->j;
246   } else {
247     p_oth = NULL;
248     pi_oth = NULL; pj_oth = NULL;
249   }
250 
251   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
252   /*-------------------------------------------------------------------*/
253   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
254   ptap->api = api;
255   api[0]    = 0;
256 
257   /* create and initialize a linked list */
258   armax    = ad->rmax+ao->rmax;
259   if (size >1) {
260     prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
261   } else {
262     prmax = p_loc->rmax;
263   }
264   nlnk_max = armax*prmax;
265   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
266   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
267 
268   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
269   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
270 
271   current_space = free_space;
272 
273   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
274   for (i=0; i<am; i++) {
275     /* diagonal portion of A */
276     nzi = adi[i+1] - adi[i];
277     for (j=0; j<nzi; j++) {
278       row  = *adj++;
279       pnz  = pi_loc[row+1] - pi_loc[row];
280       Jptr = pj_loc + pi_loc[row];
281       /* add non-zero cols of P into the sorted linked list lnk */
282       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
283     }
284     /* off-diagonal portion of A */
285     nzi = aoi[i+1] - aoi[i];
286     for (j=0; j<nzi; j++) {
287       row  = *aoj++;
288       pnz  = pi_oth[row+1] - pi_oth[row];
289       Jptr = pj_oth + pi_oth[row];
290       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
291     }
292 
293     apnz     = lnk[0];
294     api[i+1] = api[i] + apnz;
295 
296     /* if free space is not available, double the total space in the list */
297     if (current_space->local_remaining<apnz) {
298       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
299       nspacedouble++;
300     }
301 
302     /* Copy data into free space, then initialize lnk */
303     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
304     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
305 
306     current_space->array           += apnz;
307     current_space->local_used      += apnz;
308     current_space->local_remaining -= apnz;
309   }
310 
311   /* Allocate space for apj, initialize apj, and */
312   /* destroy list of free space and other temporary array(s) */
313   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
314   apj  = ptap->apj;
315   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
316   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
317 
318   /* malloc apa to store dense row A[i,:]*P */
319   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
320 
321   ptap->apa = apa;
322 
323   /* create and assemble symbolic parallel matrix Cmpi */
324   /*----------------------------------------------------*/
325   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
326   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
327   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
328 
329   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
330   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
331   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
332   for (i=0; i<am; i++) {
333     row  = i + rstart;
334     apnz = api[i+1] - api[i];
335     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
336     apj += apnz;
337   }
338   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
339   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
340 
341   ptap->destroy        = Cmpi->ops->destroy;
342   ptap->duplicate      = Cmpi->ops->duplicate;
343   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
344   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
345   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
346 
347   /* attach the supporting struct to Cmpi for reuse */
348   c       = (Mat_MPIAIJ*)Cmpi->data;
349   c->ptap = ptap;
350 
351   *C = Cmpi;
352 
353   /* set MatInfo */
354   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
355   if (afill < 1.0) afill = 1.0;
356   Cmpi->info.mallocs           = nspacedouble;
357   Cmpi->info.fill_ratio_given  = fill;
358   Cmpi->info.fill_ratio_needed = afill;
359 
360 #if defined(PETSC_USE_INFO)
361   if (api[am]) {
362     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
363     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
364   } else {
365     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
366   }
367 #endif
368   PetscFunctionReturn(0);
369 }
370 
371 #undef __FUNCT__
372 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
373 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
374 {
375   PetscErrorCode ierr;
376 
377   PetscFunctionBegin;
378   if (scall == MAT_INITIAL_MATRIX) {
379     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
380     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
381     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
382   }
383   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
384   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
385   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
386   PetscFunctionReturn(0);
387 }
388 
389 typedef struct {
390   Mat         workB;
391   PetscScalar *rvalues,*svalues;
392   MPI_Request *rwaits,*swaits;
393 } MPIAIJ_MPIDense;
394 
395 #undef __FUNCT__
396 #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
397 PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
398 {
399   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
400   PetscErrorCode  ierr;
401 
402   PetscFunctionBegin;
403   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
404   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
405   ierr = PetscFree(contents);CHKERRQ(ierr);
406   PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNCT__
410 #define __FUNCT__ "MatMatMultNumeric_MPIDense"
411 /*
412     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
413   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
414 
415   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
416 */
417 PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
418 {
419   PetscErrorCode         ierr;
420   PetscBool              flg;
421   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
422   PetscInt               nz   = aij->B->cmap->n;
423   PetscContainer         container;
424   MPIAIJ_MPIDense        *contents;
425   VecScatter             ctx   = aij->Mvctx;
426   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
427   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
428 
429   PetscFunctionBegin;
430   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
431   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
432 
433   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
434   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
435   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
436 
437   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
438 
439   ierr = PetscNew(&contents);CHKERRQ(ierr);
440   /* Create work matrix used to store off processor rows of B needed for local product */
441   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
442   /* Create work arrays needed */
443   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
444                       B->cmap->N*to->starts[to->n],&contents->svalues,
445                       from->n,&contents->rwaits,
446                       to->n,&contents->swaits);CHKERRQ(ierr);
447 
448   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
449   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
450   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
451   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
452   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
453 
454   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
455   PetscFunctionReturn(0);
456 }
457 
458 #undef __FUNCT__
459 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
460 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
461 {
462   PetscErrorCode         ierr;
463   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
464   PetscInt               nz   = aij->B->cmap->n;
465   PetscContainer         container;
466   MPIAIJ_MPIDense        *contents;
467   VecScatter             ctx   = aij->Mvctx;
468   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
469   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
470   PetscInt               m     = A->rmap->n,n=B->cmap->n;
471 
472   PetscFunctionBegin;
473   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
474   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
475   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
476   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
477   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
478   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
479   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
480 
481   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
482 
483   ierr = PetscNew(&contents);CHKERRQ(ierr);
484   /* Create work matrix used to store off processor rows of B needed for local product */
485   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
486   /* Create work arrays needed */
487   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
488                       B->cmap->N*to->starts[to->n],&contents->svalues,
489                       from->n,&contents->rwaits,
490                       to->n,&contents->swaits);CHKERRQ(ierr);
491 
492   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
493   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
494   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
495   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
496   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
497   PetscFunctionReturn(0);
498 }
499 
500 #undef __FUNCT__
501 #define __FUNCT__ "MatMPIDenseScatter"
502 /*
503     Performs an efficient scatter on the rows of B needed by this process; this is
504     a modification of the VecScatterBegin_() routines.
505 */
506 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
507 {
508   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
509   PetscErrorCode         ierr;
510   PetscScalar            *b,*w,*svalues,*rvalues;
511   VecScatter             ctx   = aij->Mvctx;
512   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
513   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
514   PetscInt               i,j,k;
515   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
516   PetscMPIInt            *sprocs,*rprocs,nrecvs;
517   MPI_Request            *swaits,*rwaits;
518   MPI_Comm               comm;
519   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
520   MPI_Status             status;
521   MPIAIJ_MPIDense        *contents;
522   PetscContainer         container;
523   Mat                    workB;
524 
525   PetscFunctionBegin;
526   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
527   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
528   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
529   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
530 
531   workB = *outworkB = contents->workB;
532   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
533   sindices = to->indices;
534   sstarts  = to->starts;
535   sprocs   = to->procs;
536   swaits   = contents->swaits;
537   svalues  = contents->svalues;
538 
539   rindices = from->indices;
540   rstarts  = from->starts;
541   rprocs   = from->procs;
542   rwaits   = contents->rwaits;
543   rvalues  = contents->rvalues;
544 
545   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
546   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
547 
548   for (i=0; i<from->n; i++) {
549     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
550   }
551 
552   for (i=0; i<to->n; i++) {
553     /* pack a message at a time */
554     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
555       for (k=0; k<ncols; k++) {
556         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
557       }
558     }
559     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
560   }
561 
562   nrecvs = from->n;
563   while (nrecvs) {
564     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
565     nrecvs--;
566     /* unpack a message at a time */
567     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
568       for (k=0; k<ncols; k++) {
569         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
570       }
571     }
572   }
573   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
574 
575   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
576   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
577   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
578   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
579   PetscFunctionReturn(0);
580 }
581 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
582 
583 #undef __FUNCT__
584 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
585 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
586 {
587   PetscErrorCode ierr;
588   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
589   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
590   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
591   Mat            workB;
592 
593   PetscFunctionBegin;
594   /* diagonal block of A times all local rows of B*/
595   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
596 
597   /* get off processor parts of B needed to complete the product */
598   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
599 
600   /* off-diagonal block of A times nonlocal rows of B */
601   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
602   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
603   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
604   PetscFunctionReturn(0);
605 }
606 
607 #undef __FUNCT__
608 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
609 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
610 {
611   PetscErrorCode ierr;
612   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
613   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
614   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
615   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
616   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
617   Mat_SeqAIJ     *p_loc,*p_oth;
618   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
619   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
620   PetscInt       cm          = C->rmap->n,anz,pnz;
621   Mat_PtAPMPI    *ptap       = c->ptap;
622   PetscScalar    *apa_sparse = ptap->apa;
623   PetscInt       *api,*apj,*apJ,i,j,k,row;
624   PetscInt       cstart = C->cmap->rstart;
625   PetscInt       cdnz,conz,k0,k1,nextp;
626   MPI_Comm       comm;
627   PetscMPIInt    size;
628 
629   PetscFunctionBegin;
630   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
631   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
632 
633   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
634   /*-----------------------------------------------------*/
635   /* update numerical values of P_oth and P_loc */
636   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
637   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
638 
639   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
640   /*----------------------------------------------------------*/
641   /* get data from symbolic products */
642   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
643   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
644   if (size >1) {
645     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
646     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
647   } else {
648     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
649   }
650 
651   api = ptap->api;
652   apj = ptap->apj;
653   for (i=0; i<cm; i++) {
654     apJ = apj + api[i];
655 
656     /* diagonal portion of A */
657     anz = adi[i+1] - adi[i];
658     adj = ad->j + adi[i];
659     ada = ad->a + adi[i];
660     for (j=0; j<anz; j++) {
661       row = adj[j];
662       pnz = pi_loc[row+1] - pi_loc[row];
663       pj  = pj_loc + pi_loc[row];
664       pa  = pa_loc + pi_loc[row];
665       /* perform sparse axpy */
666       valtmp = ada[j];
667       nextp  = 0;
668       for (k=0; nextp<pnz; k++) {
669         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
670           apa_sparse[k] += valtmp*pa[nextp++];
671         }
672       }
673       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
674     }
675 
676     /* off-diagonal portion of A */
677     anz = aoi[i+1] - aoi[i];
678     aoj = ao->j + aoi[i];
679     aoa = ao->a + aoi[i];
680     for (j=0; j<anz; j++) {
681       row = aoj[j];
682       pnz = pi_oth[row+1] - pi_oth[row];
683       pj  = pj_oth + pi_oth[row];
684       pa  = pa_oth + pi_oth[row];
685       /* perform sparse axpy */
686       valtmp = aoa[j];
687       nextp  = 0;
688       for (k=0; nextp<pnz; k++) {
689         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
690           apa_sparse[k] += valtmp*pa[nextp++];
691         }
692       }
693       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
694     }
695 
696     /* set values in C */
697     cdnz = cd->i[i+1] - cd->i[i];
698     conz = co->i[i+1] - co->i[i];
699 
700     /* 1st off-diagoanl part of C */
701     ca = coa + co->i[i];
702     k  = 0;
703     for (k0=0; k0<conz; k0++) {
704       if (apJ[k] >= cstart) break;
705       ca[k0]        = apa_sparse[k];
706       apa_sparse[k] = 0.0;
707       k++;
708     }
709 
710     /* diagonal part of C */
711     ca = cda + cd->i[i];
712     for (k1=0; k1<cdnz; k1++) {
713       ca[k1]        = apa_sparse[k];
714       apa_sparse[k] = 0.0;
715       k++;
716     }
717 
718     /* 2nd off-diagoanl part of C */
719     ca = coa + co->i[i];
720     for (; k0<conz; k0++) {
721       ca[k0]        = apa_sparse[k];
722       apa_sparse[k] = 0.0;
723       k++;
724     }
725   }
726   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
727   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
728   PetscFunctionReturn(0);
729 }
730 
731 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
732 #undef __FUNCT__
733 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
734 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
735 {
736   PetscErrorCode     ierr;
737   MPI_Comm           comm;
738   PetscMPIInt        size;
739   Mat                Cmpi;
740   Mat_PtAPMPI        *ptap;
741   PetscFreeSpaceList free_space = NULL,current_space=NULL;
742   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
743   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
744   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
745   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
746   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
747   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
748   PetscInt           nlnk_max,armax,prmax;
749   PetscReal          afill;
750   PetscScalar        *apa;
751 
752   PetscFunctionBegin;
753   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
754   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
755 
756   /* create struct Mat_PtAPMPI and attached it to C later */
757   ierr = PetscNew(&ptap);CHKERRQ(ierr);
758 
759   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
760   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
761 
762   /* get P_loc by taking all local rows of P */
763   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
764 
765   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
766   pi_loc = p_loc->i; pj_loc = p_loc->j;
767   if (size > 1) {
768     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
769     pi_oth = p_oth->i; pj_oth = p_oth->j;
770   } else {
771     p_oth  = NULL;
772     pi_oth = NULL; pj_oth = NULL;
773   }
774 
775   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
776   /*-------------------------------------------------------------------*/
777   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
778   ptap->api = api;
779   api[0]    = 0;
780 
781   /* create and initialize a linked list */
782   armax = ad->rmax+ao->rmax;
783   if (size >1) {
784     prmax = PetscMax(p_loc->rmax,p_oth->rmax);
785   } else {
786     prmax = p_loc->rmax;
787   }
788   nlnk_max = armax*prmax;
789   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
790   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
791 
792   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
793   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
794 
795   current_space = free_space;
796 
797   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
798   for (i=0; i<am; i++) {
799     /* diagonal portion of A */
800     nzi = adi[i+1] - adi[i];
801     for (j=0; j<nzi; j++) {
802       row  = *adj++;
803       pnz  = pi_loc[row+1] - pi_loc[row];
804       Jptr = pj_loc + pi_loc[row];
805       /* add non-zero cols of P into the sorted linked list lnk */
806       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
807     }
808     /* off-diagonal portion of A */
809     nzi = aoi[i+1] - aoi[i];
810     for (j=0; j<nzi; j++) {
811       row  = *aoj++;
812       pnz  = pi_oth[row+1] - pi_oth[row];
813       Jptr = pj_oth + pi_oth[row];
814       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
815     }
816 
817     apnz     = *lnk;
818     api[i+1] = api[i] + apnz;
819     if (apnz > apnz_max) apnz_max = apnz;
820 
821     /* if free space is not available, double the total space in the list */
822     if (current_space->local_remaining<apnz) {
823       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
824       nspacedouble++;
825     }
826 
827     /* Copy data into free space, then initialize lnk */
828     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
829     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
830 
831     current_space->array           += apnz;
832     current_space->local_used      += apnz;
833     current_space->local_remaining -= apnz;
834   }
835 
836   /* Allocate space for apj, initialize apj, and */
837   /* destroy list of free space and other temporary array(s) */
838   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
839   apj  = ptap->apj;
840   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
841   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
842 
843   /* create and assemble symbolic parallel matrix Cmpi */
844   /*----------------------------------------------------*/
845   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
846   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
847   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
848   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
849   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
850   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
851 
852   /* malloc apa for assembly Cmpi */
853   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
854 
855   ptap->apa = apa;
856   for (i=0; i<am; i++) {
857     row  = i + rstart;
858     apnz = api[i+1] - api[i];
859     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
860     apj += apnz;
861   }
862   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
863   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
864 
865   ptap->destroy             = Cmpi->ops->destroy;
866   ptap->duplicate           = Cmpi->ops->duplicate;
867   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
868   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
869   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
870 
871   /* attach the supporting struct to Cmpi for reuse */
872   c       = (Mat_MPIAIJ*)Cmpi->data;
873   c->ptap = ptap;
874 
875   *C = Cmpi;
876 
877   /* set MatInfo */
878   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
879   if (afill < 1.0) afill = 1.0;
880   Cmpi->info.mallocs           = nspacedouble;
881   Cmpi->info.fill_ratio_given  = fill;
882   Cmpi->info.fill_ratio_needed = afill;
883 
884 #if defined(PETSC_USE_INFO)
885   if (api[am]) {
886     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
887     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
888   } else {
889     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
890   }
891 #endif
892   PetscFunctionReturn(0);
893 }
894 
895 /*-------------------------------------------------------------------------*/
896 #undef __FUNCT__
897 #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
898 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
899 {
900   PetscErrorCode ierr;
901   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
902   PetscInt       alg=0; /* set default algorithm */
903 
904   PetscFunctionBegin;
905   if (scall == MAT_INITIAL_MATRIX) {
906     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
907     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
908     ierr = PetscOptionsEnd();CHKERRQ(ierr);
909 
910     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
911     switch (alg) {
912     case 1:
913       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
914       break;
915     case 2:
916     {
917       Mat         Pt;
918       Mat_PtAPMPI *ptap;
919       Mat_MPIAIJ  *c;
920       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
921       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
922       c        = (Mat_MPIAIJ*)(*C)->data;
923       ptap     = c->ptap;
924       ptap->Pt = Pt;
925       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
926       PetscFunctionReturn(0);
927     }
928       break;
929     default:
930       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
931       break;
932     }
933     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
934   }
935   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
936   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
937   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
938   PetscFunctionReturn(0);
939 }
940 
941 /* This routine only works when scall=MAT_REUSE_MATRIX! */
942 #undef __FUNCT__
943 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
944 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
945 {
946   PetscErrorCode ierr;
947   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
948   Mat_PtAPMPI    *ptap= c->ptap;
949   Mat            Pt=ptap->Pt;
950 
951   PetscFunctionBegin;
952   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
953   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
954   PetscFunctionReturn(0);
955 }
956 
957 /* Non-scalable version, use dense axpy */
958 #undef __FUNCT__
959 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
960 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
961 {
962   PetscErrorCode      ierr;
963   Mat_Merge_SeqsToMPI *merge;
964   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
965   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
966   Mat_PtAPMPI         *ptap;
967   PetscInt            *adj,*aJ;
968   PetscInt            i,j,k,anz,pnz,row,*cj;
969   MatScalar           *ada,*aval,*ca,valtmp;
970   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
971   MPI_Comm            comm;
972   PetscMPIInt         size,rank,taga,*len_s;
973   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
974   PetscInt            **buf_ri,**buf_rj;
975   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
976   MPI_Request         *s_waits,*r_waits;
977   MPI_Status          *status;
978   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
979   PetscInt            *ai,*aj,*coi,*coj;
980   PetscInt            *poJ,*pdJ;
981   Mat                 A_loc;
982   Mat_SeqAIJ          *a_loc;
983 
984   PetscFunctionBegin;
985   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
986   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
987   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
988 
989   ptap  = c->ptap;
990   merge = ptap->merge;
991 
992   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
993   /*--------------------------------------------------------------*/
994   /* get data from symbolic products */
995   coi  = merge->coi; coj = merge->coj;
996   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
997 
998   bi     = merge->bi; bj = merge->bj;
999   owners = merge->rowmap->range;
1000   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1001 
1002   /* get A_loc by taking all local rows of A */
1003   A_loc = ptap->A_loc;
1004   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1005   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1006   ai    = a_loc->i;
1007   aj    = a_loc->j;
1008 
1009   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
1010 
1011   for (i=0; i<am; i++) {
1012     /* 2-a) put A[i,:] to dense array aval */
1013     anz = ai[i+1] - ai[i];
1014     adj = aj + ai[i];
1015     ada = a_loc->a + ai[i];
1016     for (j=0; j<anz; j++) {
1017       aval[adj[j]] = ada[j];
1018     }
1019 
1020     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1021     /*--------------------------------------------------------------*/
1022     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1023     pnz = po->i[i+1] - po->i[i];
1024     poJ = po->j + po->i[i];
1025     pA  = po->a + po->i[i];
1026     for (j=0; j<pnz; j++) {
1027       row = poJ[j];
1028       cnz = coi[row+1] - coi[row];
1029       cj  = coj + coi[row];
1030       ca  = coa + coi[row];
1031       /* perform dense axpy */
1032       valtmp = pA[j];
1033       for (k=0; k<cnz; k++) {
1034         ca[k] += valtmp*aval[cj[k]];
1035       }
1036       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1037     }
1038 
1039     /* put the value into Cd (diagonal part) */
1040     pnz = pd->i[i+1] - pd->i[i];
1041     pdJ = pd->j + pd->i[i];
1042     pA  = pd->a + pd->i[i];
1043     for (j=0; j<pnz; j++) {
1044       row = pdJ[j];
1045       cnz = bi[row+1] - bi[row];
1046       cj  = bj + bi[row];
1047       ca  = ba + bi[row];
1048       /* perform dense axpy */
1049       valtmp = pA[j];
1050       for (k=0; k<cnz; k++) {
1051         ca[k] += valtmp*aval[cj[k]];
1052       }
1053       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1054     }
1055 
1056     /* zero the current row of Pt*A */
1057     aJ = aj + ai[i];
1058     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1059   }
1060 
1061   /* 3) send and recv matrix values coa */
1062   /*------------------------------------*/
1063   buf_ri = merge->buf_ri;
1064   buf_rj = merge->buf_rj;
1065   len_s  = merge->len_s;
1066   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1067   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1068 
1069   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1070   for (proc=0,k=0; proc<size; proc++) {
1071     if (!len_s[proc]) continue;
1072     i    = merge->owners_co[proc];
1073     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1074     k++;
1075   }
1076   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1077   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1078 
1079   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1080   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1081   ierr = PetscFree(coa);CHKERRQ(ierr);
1082 
1083   /* 4) insert local Cseq and received values into Cmpi */
1084   /*----------------------------------------------------*/
1085   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1086   for (k=0; k<merge->nrecv; k++) {
1087     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1088     nrows       = *(buf_ri_k[k]);
1089     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1090     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1091   }
1092 
1093   for (i=0; i<cm; i++) {
1094     row  = owners[rank] + i; /* global row index of C_seq */
1095     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1096     ba_i = ba + bi[i];
1097     bnz  = bi[i+1] - bi[i];
1098     /* add received vals into ba */
1099     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1100       /* i-th row */
1101       if (i == *nextrow[k]) {
1102         cnz    = *(nextci[k]+1) - *nextci[k];
1103         cj     = buf_rj[k] + *(nextci[k]);
1104         ca     = abuf_r[k] + *(nextci[k]);
1105         nextcj = 0;
1106         for (j=0; nextcj<cnz; j++) {
1107           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1108             ba_i[j] += ca[nextcj++];
1109           }
1110         }
1111         nextrow[k]++; nextci[k]++;
1112         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1113       }
1114     }
1115     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1116   }
1117   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1118   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1119 
1120   ierr = PetscFree(ba);CHKERRQ(ierr);
1121   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1122   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1123   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1124   ierr = PetscFree(aval);CHKERRQ(ierr);
1125   PetscFunctionReturn(0);
1126 }
1127 
1128 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1129 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1130 #undef __FUNCT__
1131 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
1132 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1133 {
1134   PetscErrorCode      ierr;
1135   Mat                 Cmpi,A_loc,POt,PDt;
1136   Mat_PtAPMPI         *ptap;
1137   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1138   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1139   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1140   PetscInt            nnz;
1141   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1142   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1143   PetscBT             lnkbt;
1144   MPI_Comm            comm;
1145   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1146   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1147   PetscInt            len,proc,*dnz,*onz,*owners;
1148   PetscInt            nzi,*bi,*bj;
1149   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1150   MPI_Request         *swaits,*rwaits;
1151   MPI_Status          *sstatus,rstatus;
1152   Mat_Merge_SeqsToMPI *merge;
1153   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1154   PetscReal           afill  =1.0,afill_tmp;
1155   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1156   PetscScalar         *vals;
1157   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1158 
1159   PetscFunctionBegin;
1160   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1161   /* check if matrix local sizes are compatible */
1162   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1163     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);
1164   }
1165 
1166   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1167   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1168 
1169   /* create struct Mat_PtAPMPI and attached it to C later */
1170   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1171 
1172   /* get A_loc by taking all local rows of A */
1173   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1174 
1175   ptap->A_loc = A_loc;
1176 
1177   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1178   ai    = a_loc->i;
1179   aj    = a_loc->j;
1180 
1181   /* determine symbolic Co=(p->B)^T*A - send to others */
1182   /*----------------------------------------------------*/
1183   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1184   pdt  = (Mat_SeqAIJ*)PDt->data;
1185   pdti = pdt->i; pdtj = pdt->j;
1186 
1187   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1188   pot  = (Mat_SeqAIJ*)POt->data;
1189   poti = pot->i; potj = pot->j;
1190 
1191   /* then, compute symbolic Co = (p->B)^T*A */
1192   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1193   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1194   coi[0] = 0;
1195 
1196   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1197   nnz           = fill*(poti[pon] + ai[am]);
1198   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1199   current_space = free_space;
1200 
1201   /* create and initialize a linked list */
1202   i     = PetscMax(pdt->rmax,pot->rmax);
1203   Crmax = i*a_loc->rmax*size;
1204   if (!Crmax || Crmax > aN) Crmax = aN;
1205   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1206 
1207   for (i=0; i<pon; i++) {
1208     pnz = poti[i+1] - poti[i];
1209     ptJ = potj + poti[i];
1210     for (j=0; j<pnz; j++) {
1211       row  = ptJ[j]; /* row of A_loc == col of Pot */
1212       anz  = ai[row+1] - ai[row];
1213       Jptr = aj + ai[row];
1214       /* add non-zero cols of AP into the sorted linked list lnk */
1215       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1216     }
1217     nnz = lnk[0];
1218 
1219     /* If free space is not available, double the total space in the list */
1220     if (current_space->local_remaining<nnz) {
1221       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1222       nspacedouble++;
1223     }
1224 
1225     /* Copy data into free space, and zero out denserows */
1226     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1227 
1228     current_space->array           += nnz;
1229     current_space->local_used      += nnz;
1230     current_space->local_remaining -= nnz;
1231 
1232     coi[i+1] = coi[i] + nnz;
1233   }
1234 
1235   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1236   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1237 
1238   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1239   if (afill_tmp > afill) afill = afill_tmp;
1240 
1241   /* send j-array (coj) of Co to other processors */
1242   /*----------------------------------------------*/
1243   /* determine row ownership */
1244   ierr = PetscNew(&merge);CHKERRQ(ierr);
1245   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1246 
1247   merge->rowmap->n  = pn;
1248   merge->rowmap->bs = 1;
1249 
1250   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1251   owners = merge->rowmap->range;
1252 
1253   /* determine the number of messages to send, their lengths */
1254   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1255   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
1256 
1257   len_s        = merge->len_s;
1258   merge->nsend = 0;
1259 
1260   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1261   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1262 
1263   proc = 0;
1264   for (i=0; i<pon; i++) {
1265     while (prmap[i] >= owners[proc+1]) proc++;
1266     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1267     len_s[proc] += coi[i+1] - coi[i];
1268   }
1269 
1270   len          = 0; /* max length of buf_si[] */
1271   owners_co[0] = 0;
1272   for (proc=0; proc<size; proc++) {
1273     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1274     if (len_si[proc]) {
1275       merge->nsend++;
1276       len_si[proc] = 2*(len_si[proc] + 1);
1277       len         += len_si[proc];
1278     }
1279   }
1280 
1281   /* determine the number and length of messages to receive for coi and coj  */
1282   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1283   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1284 
1285   /* post the Irecv and Isend of coj */
1286   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1287   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1288   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1289   for (proc=0, k=0; proc<size; proc++) {
1290     if (!len_s[proc]) continue;
1291     i    = owners_co[proc];
1292     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1293     k++;
1294   }
1295 
1296   /* receives and sends of coj are complete */
1297   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1298   for (i=0; i<merge->nrecv; i++) {
1299     PetscMPIInt icompleted;
1300     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1301   }
1302   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1303   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1304 
1305   /* send and recv coi */
1306   /*-------------------*/
1307   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1308   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1309   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1310   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1311   for (proc=0,k=0; proc<size; proc++) {
1312     if (!len_s[proc]) continue;
1313     /* form outgoing message for i-structure:
1314          buf_si[0]:                 nrows to be sent
1315                [1:nrows]:           row index (global)
1316                [nrows+1:2*nrows+1]: i-structure index
1317     */
1318     /*-------------------------------------------*/
1319     nrows       = len_si[proc]/2 - 1;
1320     buf_si_i    = buf_si + nrows+1;
1321     buf_si[0]   = nrows;
1322     buf_si_i[0] = 0;
1323     nrows       = 0;
1324     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1325       nzi               = coi[i+1] - coi[i];
1326       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1327       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1328       nrows++;
1329     }
1330     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1331     k++;
1332     buf_si += len_si[proc];
1333   }
1334   i = merge->nrecv;
1335   while (i--) {
1336     PetscMPIInt icompleted;
1337     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1338   }
1339   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1340   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1341   ierr = PetscFree(len_si);CHKERRQ(ierr);
1342   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1343   ierr = PetscFree(swaits);CHKERRQ(ierr);
1344   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1345   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1346 
1347   /* compute the local portion of C (mpi mat) */
1348   /*------------------------------------------*/
1349   /* allocate bi array and free space for accumulating nonzero column info */
1350   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1351   bi[0] = 0;
1352 
1353   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1354   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1355   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1356   current_space = free_space;
1357 
1358   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1359   for (k=0; k<merge->nrecv; k++) {
1360     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1361     nrows       = *buf_ri_k[k];
1362     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1363     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1364   }
1365 
1366   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1367   rmax = 0;
1368   for (i=0; i<pn; i++) {
1369     /* add pdt[i,:]*AP into lnk */
1370     pnz = pdti[i+1] - pdti[i];
1371     ptJ = pdtj + pdti[i];
1372     for (j=0; j<pnz; j++) {
1373       row  = ptJ[j];  /* row of AP == col of Pt */
1374       anz  = ai[row+1] - ai[row];
1375       Jptr = aj + ai[row];
1376       /* add non-zero cols of AP into the sorted linked list lnk */
1377       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1378     }
1379 
1380     /* add received col data into lnk */
1381     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1382       if (i == *nextrow[k]) { /* i-th row */
1383         nzi  = *(nextci[k]+1) - *nextci[k];
1384         Jptr = buf_rj[k] + *nextci[k];
1385         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1386         nextrow[k]++; nextci[k]++;
1387       }
1388     }
1389     nnz = lnk[0];
1390 
1391     /* if free space is not available, make more free space */
1392     if (current_space->local_remaining<nnz) {
1393       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1394       nspacedouble++;
1395     }
1396     /* copy data into free space, then initialize lnk */
1397     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1398     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1399 
1400     current_space->array           += nnz;
1401     current_space->local_used      += nnz;
1402     current_space->local_remaining -= nnz;
1403 
1404     bi[i+1] = bi[i] + nnz;
1405     if (nnz > rmax) rmax = nnz;
1406   }
1407   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1408 
1409   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1410   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1411 
1412   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1413   if (afill_tmp > afill) afill = afill_tmp;
1414   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
1415   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1416   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1417 
1418   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1419   /*----------------------------------------------------------------------------------*/
1420   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1421 
1422   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1423   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1424   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1425   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1426   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1427   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1428   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1429   for (i=0; i<pn; i++) {
1430     row  = i + rstart;
1431     nnz  = bi[i+1] - bi[i];
1432     Jptr = bj + bi[i];
1433     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1434   }
1435   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1436   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1437   ierr = PetscFree(vals);CHKERRQ(ierr);
1438 
1439   merge->bi        = bi;
1440   merge->bj        = bj;
1441   merge->coi       = coi;
1442   merge->coj       = coj;
1443   merge->buf_ri    = buf_ri;
1444   merge->buf_rj    = buf_rj;
1445   merge->owners_co = owners_co;
1446   merge->destroy   = Cmpi->ops->destroy;
1447   merge->duplicate = Cmpi->ops->duplicate;
1448 
1449   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1450   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1451   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1452 
1453   /* attach the supporting struct to Cmpi for reuse */
1454   c           = (Mat_MPIAIJ*)Cmpi->data;
1455   c->ptap     = ptap;
1456   ptap->api   = NULL;
1457   ptap->apj   = NULL;
1458   ptap->merge = merge;
1459   ptap->rmax  = rmax;
1460 
1461   *C = Cmpi;
1462 #if defined(PETSC_USE_INFO)
1463   if (bi[pn] != 0) {
1464     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
1465     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1466   } else {
1467     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1468   }
1469 #endif
1470   PetscFunctionReturn(0);
1471 }
1472 
1473 #undef __FUNCT__
1474 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
1475 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1476 {
1477   PetscErrorCode      ierr;
1478   Mat_Merge_SeqsToMPI *merge;
1479   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1480   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1481   Mat_PtAPMPI         *ptap;
1482   PetscInt            *adj;
1483   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1484   MatScalar           *ada,*ca,valtmp;
1485   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1486   MPI_Comm            comm;
1487   PetscMPIInt         size,rank,taga,*len_s;
1488   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1489   PetscInt            **buf_ri,**buf_rj;
1490   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1491   MPI_Request         *s_waits,*r_waits;
1492   MPI_Status          *status;
1493   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1494   PetscInt            *ai,*aj,*coi,*coj;
1495   PetscInt            *poJ,*pdJ;
1496   Mat                 A_loc;
1497   Mat_SeqAIJ          *a_loc;
1498 
1499   PetscFunctionBegin;
1500   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1501   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1502   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1503 
1504   ptap  = c->ptap;
1505   merge = ptap->merge;
1506 
1507   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1508   /*------------------------------------------*/
1509   /* get data from symbolic products */
1510   coi    = merge->coi; coj = merge->coj;
1511   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1512   bi     = merge->bi; bj = merge->bj;
1513   owners = merge->rowmap->range;
1514   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1515 
1516   /* get A_loc by taking all local rows of A */
1517   A_loc = ptap->A_loc;
1518   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1519   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1520   ai    = a_loc->i;
1521   aj    = a_loc->j;
1522 
1523   for (i=0; i<am; i++) {
1524     anz = ai[i+1] - ai[i];
1525     adj = aj + ai[i];
1526     ada = a_loc->a + ai[i];
1527 
1528     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1529     /*-------------------------------------------------------------*/
1530     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1531     pnz = po->i[i+1] - po->i[i];
1532     poJ = po->j + po->i[i];
1533     pA  = po->a + po->i[i];
1534     for (j=0; j<pnz; j++) {
1535       row = poJ[j];
1536       cj  = coj + coi[row];
1537       ca  = coa + coi[row];
1538       /* perform sparse axpy */
1539       nexta  = 0;
1540       valtmp = pA[j];
1541       for (k=0; nexta<anz; k++) {
1542         if (cj[k] == adj[nexta]) {
1543           ca[k] += valtmp*ada[nexta];
1544           nexta++;
1545         }
1546       }
1547       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1548     }
1549 
1550     /* put the value into Cd (diagonal part) */
1551     pnz = pd->i[i+1] - pd->i[i];
1552     pdJ = pd->j + pd->i[i];
1553     pA  = pd->a + pd->i[i];
1554     for (j=0; j<pnz; j++) {
1555       row = pdJ[j];
1556       cj  = bj + bi[row];
1557       ca  = ba + bi[row];
1558       /* perform sparse axpy */
1559       nexta  = 0;
1560       valtmp = pA[j];
1561       for (k=0; nexta<anz; k++) {
1562         if (cj[k] == adj[nexta]) {
1563           ca[k] += valtmp*ada[nexta];
1564           nexta++;
1565         }
1566       }
1567       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1568     }
1569   }
1570 
1571   /* 3) send and recv matrix values coa */
1572   /*------------------------------------*/
1573   buf_ri = merge->buf_ri;
1574   buf_rj = merge->buf_rj;
1575   len_s  = merge->len_s;
1576   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1577   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1578 
1579   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1580   for (proc=0,k=0; proc<size; proc++) {
1581     if (!len_s[proc]) continue;
1582     i    = merge->owners_co[proc];
1583     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1584     k++;
1585   }
1586   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1587   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1588 
1589   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1590   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1591   ierr = PetscFree(coa);CHKERRQ(ierr);
1592 
1593   /* 4) insert local Cseq and received values into Cmpi */
1594   /*----------------------------------------------------*/
1595   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1596   for (k=0; k<merge->nrecv; k++) {
1597     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1598     nrows       = *(buf_ri_k[k]);
1599     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1600     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1601   }
1602 
1603   for (i=0; i<cm; i++) {
1604     row  = owners[rank] + i; /* global row index of C_seq */
1605     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1606     ba_i = ba + bi[i];
1607     bnz  = bi[i+1] - bi[i];
1608     /* add received vals into ba */
1609     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1610       /* i-th row */
1611       if (i == *nextrow[k]) {
1612         cnz    = *(nextci[k]+1) - *nextci[k];
1613         cj     = buf_rj[k] + *(nextci[k]);
1614         ca     = abuf_r[k] + *(nextci[k]);
1615         nextcj = 0;
1616         for (j=0; nextcj<cnz; j++) {
1617           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1618             ba_i[j] += ca[nextcj++];
1619           }
1620         }
1621         nextrow[k]++; nextci[k]++;
1622         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1623       }
1624     }
1625     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1626   }
1627   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1628   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1629 
1630   ierr = PetscFree(ba);CHKERRQ(ierr);
1631   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1632   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1633   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1634   PetscFunctionReturn(0);
1635 }
1636 
1637 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1638 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
1639    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1640 #undef __FUNCT__
1641 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
1642 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1643 {
1644   PetscErrorCode      ierr;
1645   Mat                 Cmpi,A_loc,POt,PDt;
1646   Mat_PtAPMPI         *ptap;
1647   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1648   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1649   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1650   PetscInt            nnz;
1651   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1652   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1653   MPI_Comm            comm;
1654   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1655   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1656   PetscInt            len,proc,*dnz,*onz,*owners;
1657   PetscInt            nzi,*bi,*bj;
1658   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1659   MPI_Request         *swaits,*rwaits;
1660   MPI_Status          *sstatus,rstatus;
1661   Mat_Merge_SeqsToMPI *merge;
1662   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1663   PetscReal           afill  =1.0,afill_tmp;
1664   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1665   PetscScalar         *vals;
1666   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1667 
1668   PetscFunctionBegin;
1669   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1670   /* check if matrix local sizes are compatible */
1671   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1672     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);
1673   }
1674 
1675   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1676   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1677 
1678   /* create struct Mat_PtAPMPI and attached it to C later */
1679   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1680 
1681   /* get A_loc by taking all local rows of A */
1682   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1683 
1684   ptap->A_loc = A_loc;
1685   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1686   ai          = a_loc->i;
1687   aj          = a_loc->j;
1688 
1689   /* determine symbolic Co=(p->B)^T*A - send to others */
1690   /*----------------------------------------------------*/
1691   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1692   pdt  = (Mat_SeqAIJ*)PDt->data;
1693   pdti = pdt->i; pdtj = pdt->j;
1694 
1695   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1696   pot  = (Mat_SeqAIJ*)POt->data;
1697   poti = pot->i; potj = pot->j;
1698 
1699   /* then, compute symbolic Co = (p->B)^T*A */
1700   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1701                          >= (num of nonzero rows of C_seq) - pn */
1702   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1703   coi[0] = 0;
1704 
1705   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1706   nnz           = fill*(poti[pon] + ai[am]);
1707   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1708   current_space = free_space;
1709 
1710   /* create and initialize a linked list */
1711   i     = PetscMax(pdt->rmax,pot->rmax);
1712   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1713   if (!Crmax || Crmax > aN) Crmax = aN;
1714   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1715 
1716   for (i=0; i<pon; i++) {
1717     pnz = poti[i+1] - poti[i];
1718     ptJ = potj + poti[i];
1719     for (j=0; j<pnz; j++) {
1720       row  = ptJ[j]; /* row of A_loc == col of Pot */
1721       anz  = ai[row+1] - ai[row];
1722       Jptr = aj + ai[row];
1723       /* add non-zero cols of AP into the sorted linked list lnk */
1724       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1725     }
1726     nnz = lnk[0];
1727 
1728     /* If free space is not available, double the total space in the list */
1729     if (current_space->local_remaining<nnz) {
1730       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1731       nspacedouble++;
1732     }
1733 
1734     /* Copy data into free space, and zero out denserows */
1735     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1736 
1737     current_space->array           += nnz;
1738     current_space->local_used      += nnz;
1739     current_space->local_remaining -= nnz;
1740 
1741     coi[i+1] = coi[i] + nnz;
1742   }
1743 
1744   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1745   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1746 
1747   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1748   if (afill_tmp > afill) afill = afill_tmp;
1749 
1750   /* send j-array (coj) of Co to other processors */
1751   /*----------------------------------------------*/
1752   /* determine row ownership */
1753   ierr = PetscNew(&merge);CHKERRQ(ierr);
1754   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1755 
1756   merge->rowmap->n  = pn;
1757   merge->rowmap->bs = 1;
1758 
1759   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1760   owners = merge->rowmap->range;
1761 
1762   /* determine the number of messages to send, their lengths */
1763   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1764   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
1765 
1766   len_s        = merge->len_s;
1767   merge->nsend = 0;
1768 
1769   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1770   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1771 
1772   proc = 0;
1773   for (i=0; i<pon; i++) {
1774     while (prmap[i] >= owners[proc+1]) proc++;
1775     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1776     len_s[proc] += coi[i+1] - coi[i];
1777   }
1778 
1779   len          = 0; /* max length of buf_si[] */
1780   owners_co[0] = 0;
1781   for (proc=0; proc<size; proc++) {
1782     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1783     if (len_si[proc]) {
1784       merge->nsend++;
1785       len_si[proc] = 2*(len_si[proc] + 1);
1786       len         += len_si[proc];
1787     }
1788   }
1789 
1790   /* determine the number and length of messages to receive for coi and coj  */
1791   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1792   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1793 
1794   /* post the Irecv and Isend of coj */
1795   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1796   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1797   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1798   for (proc=0, k=0; proc<size; proc++) {
1799     if (!len_s[proc]) continue;
1800     i    = owners_co[proc];
1801     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1802     k++;
1803   }
1804 
1805   /* receives and sends of coj are complete */
1806   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1807   for (i=0; i<merge->nrecv; i++) {
1808     PetscMPIInt icompleted;
1809     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1810   }
1811   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1812   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1813 
1814   /* send and recv coi */
1815   /*-------------------*/
1816   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1817   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1818   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1819   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1820   for (proc=0,k=0; proc<size; proc++) {
1821     if (!len_s[proc]) continue;
1822     /* form outgoing message for i-structure:
1823          buf_si[0]:                 nrows to be sent
1824                [1:nrows]:           row index (global)
1825                [nrows+1:2*nrows+1]: i-structure index
1826     */
1827     /*-------------------------------------------*/
1828     nrows       = len_si[proc]/2 - 1;
1829     buf_si_i    = buf_si + nrows+1;
1830     buf_si[0]   = nrows;
1831     buf_si_i[0] = 0;
1832     nrows       = 0;
1833     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1834       nzi               = coi[i+1] - coi[i];
1835       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1836       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1837       nrows++;
1838     }
1839     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1840     k++;
1841     buf_si += len_si[proc];
1842   }
1843   i = merge->nrecv;
1844   while (i--) {
1845     PetscMPIInt icompleted;
1846     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1847   }
1848   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1849   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1850   ierr = PetscFree(len_si);CHKERRQ(ierr);
1851   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1852   ierr = PetscFree(swaits);CHKERRQ(ierr);
1853   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1854   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1855 
1856   /* compute the local portion of C (mpi mat) */
1857   /*------------------------------------------*/
1858   /* allocate bi array and free space for accumulating nonzero column info */
1859   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1860   bi[0] = 0;
1861 
1862   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1863   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1864   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1865   current_space = free_space;
1866 
1867   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1868   for (k=0; k<merge->nrecv; k++) {
1869     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1870     nrows       = *buf_ri_k[k];
1871     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1872     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1873   }
1874 
1875   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1876   rmax = 0;
1877   for (i=0; i<pn; i++) {
1878     /* add pdt[i,:]*AP into lnk */
1879     pnz = pdti[i+1] - pdti[i];
1880     ptJ = pdtj + pdti[i];
1881     for (j=0; j<pnz; j++) {
1882       row  = ptJ[j];  /* row of AP == col of Pt */
1883       anz  = ai[row+1] - ai[row];
1884       Jptr = aj + ai[row];
1885       /* add non-zero cols of AP into the sorted linked list lnk */
1886       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1887     }
1888 
1889     /* add received col data into lnk */
1890     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1891       if (i == *nextrow[k]) { /* i-th row */
1892         nzi  = *(nextci[k]+1) - *nextci[k];
1893         Jptr = buf_rj[k] + *nextci[k];
1894         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1895         nextrow[k]++; nextci[k]++;
1896       }
1897     }
1898     nnz = lnk[0];
1899 
1900     /* if free space is not available, make more free space */
1901     if (current_space->local_remaining<nnz) {
1902       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1903       nspacedouble++;
1904     }
1905     /* copy data into free space, then initialize lnk */
1906     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1907     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1908 
1909     current_space->array           += nnz;
1910     current_space->local_used      += nnz;
1911     current_space->local_remaining -= nnz;
1912 
1913     bi[i+1] = bi[i] + nnz;
1914     if (nnz > rmax) rmax = nnz;
1915   }
1916   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1917 
1918   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1919   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1920   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1921   if (afill_tmp > afill) afill = afill_tmp;
1922   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1923   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1924   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1925 
1926   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1927   /*----------------------------------------------------------------------------------*/
1928   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1929 
1930   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1931   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1932   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1933   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1934   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1935   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1936   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1937   for (i=0; i<pn; i++) {
1938     row  = i + rstart;
1939     nnz  = bi[i+1] - bi[i];
1940     Jptr = bj + bi[i];
1941     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1942   }
1943   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1944   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1945   ierr = PetscFree(vals);CHKERRQ(ierr);
1946 
1947   merge->bi        = bi;
1948   merge->bj        = bj;
1949   merge->coi       = coi;
1950   merge->coj       = coj;
1951   merge->buf_ri    = buf_ri;
1952   merge->buf_rj    = buf_rj;
1953   merge->owners_co = owners_co;
1954   merge->destroy   = Cmpi->ops->destroy;
1955   merge->duplicate = Cmpi->ops->duplicate;
1956 
1957   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1958   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1959   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1960 
1961   /* attach the supporting struct to Cmpi for reuse */
1962   c = (Mat_MPIAIJ*)Cmpi->data;
1963 
1964   c->ptap     = ptap;
1965   ptap->api   = NULL;
1966   ptap->apj   = NULL;
1967   ptap->merge = merge;
1968   ptap->rmax  = rmax;
1969   ptap->apa   = NULL;
1970 
1971   *C = Cmpi;
1972 #if defined(PETSC_USE_INFO)
1973   if (bi[pn] != 0) {
1974     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
1975     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1976   } else {
1977     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1978   }
1979 #endif
1980   PetscFunctionReturn(0);
1981 }
1982