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