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