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