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