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