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