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