xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 42ce410b93811b6a1b27796b46869d3e6c55267a)
1 /*
2    Basic functions for basic parallel dense matrices.
3    Portions of this code are under:
4    Copyright (c) 2022 Advanced Micro Devices, Inc. All rights reserved.
5 */
6 
7 #include <../src/mat/impls/dense/mpi/mpidense.h> /*I   "petscmat.h"  I*/
8 #include <../src/mat/impls/aij/mpi/mpiaij.h>
9 #include <petscblaslapack.h>
10 #include <petsc/private/vecimpl.h>
11 
12 /*@
13   MatDenseGetLocalMatrix - For a `MATMPIDENSE` or `MATSEQDENSE` matrix returns the sequential
14   matrix that represents the operator. For sequential matrices it returns itself.
15 
16   Input Parameter:
17 . A - the sequential or MPI `MATDENSE` matrix
18 
19   Output Parameter:
20 . B - the inner matrix
21 
22   Level: intermediate
23 
24 .seealso: [](ch_matrices), `Mat`, `MATDENSE`, `MATMPIDENSE`, `MATSEQDENSE`
25 @*/
26 PetscErrorCode MatDenseGetLocalMatrix(Mat A, Mat *B)
27 {
28   Mat_MPIDense *mat = (Mat_MPIDense *)A->data;
29   PetscBool     flg;
30 
31   PetscFunctionBegin;
32   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
33   PetscAssertPointer(B, 2);
34   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIDENSE, &flg));
35   if (flg) *B = mat->A;
36   else {
37     PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQDENSE, &flg));
38     PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for matrix type %s", ((PetscObject)A)->type_name);
39     *B = A;
40   }
41   PetscFunctionReturn(PETSC_SUCCESS);
42 }
43 
44 static PetscErrorCode MatCopy_MPIDense(Mat A, Mat B, MatStructure s)
45 {
46   Mat_MPIDense *Amat = (Mat_MPIDense *)A->data;
47   Mat_MPIDense *Bmat = (Mat_MPIDense *)B->data;
48 
49   PetscFunctionBegin;
50   PetscCall(MatCopy(Amat->A, Bmat->A, s));
51   PetscFunctionReturn(PETSC_SUCCESS);
52 }
53 
54 PetscErrorCode MatShift_MPIDense(Mat A, PetscScalar alpha)
55 {
56   Mat_MPIDense *mat = (Mat_MPIDense *)A->data;
57   PetscInt      j, lda, rstart = A->rmap->rstart, rend = A->rmap->rend, rend2;
58   PetscScalar  *v;
59 
60   PetscFunctionBegin;
61   PetscCall(MatDenseGetArray(mat->A, &v));
62   PetscCall(MatDenseGetLDA(mat->A, &lda));
63   rend2 = PetscMin(rend, A->cmap->N);
64   if (rend2 > rstart) {
65     for (j = rstart; j < rend2; j++) v[j - rstart + j * lda] += alpha;
66     PetscCall(PetscLogFlops(rend2 - rstart));
67   }
68   PetscCall(MatDenseRestoreArray(mat->A, &v));
69   PetscFunctionReturn(PETSC_SUCCESS);
70 }
71 
72 static PetscErrorCode MatGetRow_MPIDense(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
73 {
74   Mat_MPIDense *mat = (Mat_MPIDense *)A->data;
75   PetscInt      lrow, rstart = A->rmap->rstart, rend = A->rmap->rend;
76 
77   PetscFunctionBegin;
78   PetscCheck(row >= rstart && row < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "only local rows");
79   lrow = row - rstart;
80   PetscCall(MatGetRow(mat->A, lrow, nz, (const PetscInt **)idx, (const PetscScalar **)v));
81   PetscFunctionReturn(PETSC_SUCCESS);
82 }
83 
84 static PetscErrorCode MatRestoreRow_MPIDense(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
85 {
86   Mat_MPIDense *mat = (Mat_MPIDense *)A->data;
87   PetscInt      lrow, rstart = A->rmap->rstart, rend = A->rmap->rend;
88 
89   PetscFunctionBegin;
90   PetscCheck(row >= rstart && row < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "only local rows");
91   lrow = row - rstart;
92   PetscCall(MatRestoreRow(mat->A, lrow, nz, (const PetscInt **)idx, (const PetscScalar **)v));
93   PetscFunctionReturn(PETSC_SUCCESS);
94 }
95 
96 static PetscErrorCode MatGetDiagonalBlock_MPIDense(Mat A, Mat *a)
97 {
98   Mat_MPIDense *mdn = (Mat_MPIDense *)A->data;
99   PetscInt      m = A->rmap->n, rstart = A->rmap->rstart;
100   PetscScalar  *array;
101   MPI_Comm      comm;
102   PetscBool     flg;
103   Mat           B;
104 
105   PetscFunctionBegin;
106   PetscCall(MatHasCongruentLayouts(A, &flg));
107   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only square matrices supported.");
108   PetscCall(PetscObjectQuery((PetscObject)A, "DiagonalBlock", (PetscObject *)&B));
109   if (!B) { /* This should use MatDenseGetSubMatrix (not create), but we would need a call like MatRestoreDiagonalBlock */
110 #if PetscDefined(HAVE_CUDA)
111     PetscCall(PetscObjectTypeCompare((PetscObject)mdn->A, MATSEQDENSECUDA, &flg));
112     PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not coded for %s. Send an email to petsc-dev@mcs.anl.gov to request this feature", MATSEQDENSECUDA);
113 #elif PetscDefined(HAVE_HIP)
114     PetscCall(PetscObjectTypeCompare((PetscObject)mdn->A, MATSEQDENSEHIP, &flg));
115     PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not coded for %s. Send an email to petsc-dev@mcs.anl.gov to request this feature", MATSEQDENSEHIP);
116 #endif
117     PetscCall(PetscObjectGetComm((PetscObject)mdn->A, &comm));
118     PetscCall(MatCreate(comm, &B));
119     PetscCall(MatSetSizes(B, m, m, m, m));
120     PetscCall(MatSetType(B, ((PetscObject)mdn->A)->type_name));
121     PetscCall(MatDenseGetArrayRead(mdn->A, (const PetscScalar **)&array));
122     PetscCall(MatSeqDenseSetPreallocation(B, array + m * rstart));
123     PetscCall(MatDenseRestoreArrayRead(mdn->A, (const PetscScalar **)&array));
124     PetscCall(PetscObjectCompose((PetscObject)A, "DiagonalBlock", (PetscObject)B));
125     *a = B;
126     PetscCall(MatDestroy(&B));
127   } else *a = B;
128   PetscFunctionReturn(PETSC_SUCCESS);
129 }
130 
131 static PetscErrorCode MatSetValues_MPIDense(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], const PetscScalar v[], InsertMode addv)
132 {
133   Mat_MPIDense *A = (Mat_MPIDense *)mat->data;
134   PetscInt      i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend, row;
135   PetscBool     roworiented = A->roworiented;
136 
137   PetscFunctionBegin;
138   for (i = 0; i < m; i++) {
139     if (idxm[i] < 0) continue;
140     PetscCheck(idxm[i] < mat->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large");
141     if (idxm[i] >= rstart && idxm[i] < rend) {
142       row = idxm[i] - rstart;
143       if (roworiented) {
144         PetscCall(MatSetValues(A->A, 1, &row, n, idxn, PetscSafePointerPlusOffset(v, i * n), addv));
145       } else {
146         for (j = 0; j < n; j++) {
147           if (idxn[j] < 0) continue;
148           PetscCheck(idxn[j] < mat->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large");
149           PetscCall(MatSetValues(A->A, 1, &row, 1, &idxn[j], PetscSafePointerPlusOffset(v, i + j * m), addv));
150         }
151       }
152     } else if (!A->donotstash) {
153       mat->assembled = PETSC_FALSE;
154       if (roworiented) {
155         PetscCall(MatStashValuesRow_Private(&mat->stash, idxm[i], n, idxn, PetscSafePointerPlusOffset(v, i * n), PETSC_FALSE));
156       } else {
157         PetscCall(MatStashValuesCol_Private(&mat->stash, idxm[i], n, idxn, PetscSafePointerPlusOffset(v, i), m, PETSC_FALSE));
158       }
159     }
160   }
161   PetscFunctionReturn(PETSC_SUCCESS);
162 }
163 
164 static PetscErrorCode MatGetValues_MPIDense(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[])
165 {
166   Mat_MPIDense *mdn = (Mat_MPIDense *)mat->data;
167   PetscInt      i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend, row;
168 
169   PetscFunctionBegin;
170   for (i = 0; i < m; i++) {
171     if (idxm[i] < 0) continue; /* negative row */
172     PetscCheck(idxm[i] < mat->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large");
173     if (idxm[i] >= rstart && idxm[i] < rend) {
174       row = idxm[i] - rstart;
175       for (j = 0; j < n; j++) {
176         if (idxn[j] < 0) continue; /* negative column */
177         PetscCheck(idxn[j] < mat->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large");
178         PetscCall(MatGetValues(mdn->A, 1, &row, 1, &idxn[j], v + i * n + j));
179       }
180     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only local values currently supported");
181   }
182   PetscFunctionReturn(PETSC_SUCCESS);
183 }
184 
185 static PetscErrorCode MatDenseGetLDA_MPIDense(Mat A, PetscInt *lda)
186 {
187   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
188 
189   PetscFunctionBegin;
190   PetscCall(MatDenseGetLDA(a->A, lda));
191   PetscFunctionReturn(PETSC_SUCCESS);
192 }
193 
194 static PetscErrorCode MatDenseSetLDA_MPIDense(Mat A, PetscInt lda)
195 {
196   Mat_MPIDense *a     = (Mat_MPIDense *)A->data;
197   MatType       mtype = MATSEQDENSE;
198 
199   PetscFunctionBegin;
200   if (!a->A) {
201     PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
202     PetscCall(PetscLayoutSetUp(A->rmap));
203     PetscCall(PetscLayoutSetUp(A->cmap));
204     PetscCall(MatCreate(PETSC_COMM_SELF, &a->A));
205     PetscCall(MatSetSizes(a->A, A->rmap->n, A->cmap->N, A->rmap->n, A->cmap->N));
206 #if PetscDefined(HAVE_CUDA)
207     PetscBool iscuda;
208     PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIDENSECUDA, &iscuda));
209     if (iscuda) mtype = MATSEQDENSECUDA;
210 #elif PetscDefined(HAVE_HIP)
211     PetscBool iship;
212     PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIDENSEHIP, &iship));
213     if (iship) mtype = MATSEQDENSEHIP;
214 #endif
215     PetscCall(MatSetType(a->A, mtype));
216   }
217   PetscCall(MatDenseSetLDA(a->A, lda));
218   PetscFunctionReturn(PETSC_SUCCESS);
219 }
220 
221 static PetscErrorCode MatDenseGetArray_MPIDense(Mat A, PetscScalar **array)
222 {
223   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
224 
225   PetscFunctionBegin;
226   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
227   PetscCall(MatDenseGetArray(a->A, array));
228   PetscFunctionReturn(PETSC_SUCCESS);
229 }
230 
231 static PetscErrorCode MatDenseGetArrayRead_MPIDense(Mat A, PetscScalar **array)
232 {
233   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
234 
235   PetscFunctionBegin;
236   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
237   PetscCall(MatDenseGetArrayRead(a->A, (const PetscScalar **)array));
238   PetscFunctionReturn(PETSC_SUCCESS);
239 }
240 
241 static PetscErrorCode MatDenseGetArrayWrite_MPIDense(Mat A, PetscScalar **array)
242 {
243   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
244 
245   PetscFunctionBegin;
246   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
247   PetscCall(MatDenseGetArrayWrite(a->A, array));
248   PetscFunctionReturn(PETSC_SUCCESS);
249 }
250 
251 static PetscErrorCode MatDensePlaceArray_MPIDense(Mat A, const PetscScalar *array)
252 {
253   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
254 
255   PetscFunctionBegin;
256   PetscCheck(!a->vecinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first");
257   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
258   PetscCall(MatDensePlaceArray(a->A, array));
259   PetscFunctionReturn(PETSC_SUCCESS);
260 }
261 
262 static PetscErrorCode MatDenseResetArray_MPIDense(Mat A)
263 {
264   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
265 
266   PetscFunctionBegin;
267   PetscCheck(!a->vecinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first");
268   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
269   PetscCall(MatDenseResetArray(a->A));
270   PetscFunctionReturn(PETSC_SUCCESS);
271 }
272 
273 static PetscErrorCode MatDenseReplaceArray_MPIDense(Mat A, const PetscScalar *array)
274 {
275   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
276 
277   PetscFunctionBegin;
278   PetscCheck(!a->vecinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first");
279   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
280   PetscCall(MatDenseReplaceArray(a->A, array));
281   PetscFunctionReturn(PETSC_SUCCESS);
282 }
283 
284 static PetscErrorCode MatCreateSubMatrix_MPIDense(Mat A, IS isrow, IS iscol, MatReuse scall, Mat *B)
285 {
286   Mat_MPIDense      *mat = (Mat_MPIDense *)A->data, *newmatd;
287   PetscInt           lda, i, j, rstart, rend, nrows, ncols, Ncols, nlrows, nlcols;
288   const PetscInt    *irow, *icol;
289   const PetscScalar *v;
290   PetscScalar       *bv;
291   Mat                newmat;
292   IS                 iscol_local;
293   MPI_Comm           comm_is, comm_mat;
294 
295   PetscFunctionBegin;
296   PetscCall(PetscObjectGetComm((PetscObject)A, &comm_mat));
297   PetscCall(PetscObjectGetComm((PetscObject)iscol, &comm_is));
298   PetscCheck(comm_mat == comm_is, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMECOMM, "IS communicator must match matrix communicator");
299 
300   PetscCall(ISAllGather(iscol, &iscol_local));
301   PetscCall(ISGetIndices(isrow, &irow));
302   PetscCall(ISGetIndices(iscol_local, &icol));
303   PetscCall(ISGetLocalSize(isrow, &nrows));
304   PetscCall(ISGetLocalSize(iscol, &ncols));
305   PetscCall(ISGetSize(iscol, &Ncols)); /* global number of columns, size of iscol_local */
306 
307   /* No parallel redistribution currently supported! Should really check each index set
308      to confirm that it is OK.  ... Currently supports only submatrix same partitioning as
309      original matrix! */
310 
311   PetscCall(MatGetLocalSize(A, &nlrows, &nlcols));
312   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
313 
314   /* Check submatrix call */
315   if (scall == MAT_REUSE_MATRIX) {
316     /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
317     /* Really need to test rows and column sizes! */
318     newmat = *B;
319   } else {
320     /* Create and fill new matrix */
321     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &newmat));
322     PetscCall(MatSetSizes(newmat, nrows, ncols, PETSC_DECIDE, Ncols));
323     PetscCall(MatSetType(newmat, ((PetscObject)A)->type_name));
324     PetscCall(MatMPIDenseSetPreallocation(newmat, NULL));
325   }
326 
327   /* Now extract the data pointers and do the copy, column at a time */
328   newmatd = (Mat_MPIDense *)newmat->data;
329   PetscCall(MatDenseGetArray(newmatd->A, &bv));
330   PetscCall(MatDenseGetArrayRead(mat->A, &v));
331   PetscCall(MatDenseGetLDA(mat->A, &lda));
332   for (i = 0; i < Ncols; i++) {
333     const PetscScalar *av = v + lda * icol[i];
334     for (j = 0; j < nrows; j++) *bv++ = av[irow[j] - rstart];
335   }
336   PetscCall(MatDenseRestoreArrayRead(mat->A, &v));
337   PetscCall(MatDenseRestoreArray(newmatd->A, &bv));
338 
339   /* Assemble the matrices so that the correct flags are set */
340   PetscCall(MatAssemblyBegin(newmat, MAT_FINAL_ASSEMBLY));
341   PetscCall(MatAssemblyEnd(newmat, MAT_FINAL_ASSEMBLY));
342 
343   /* Free work space */
344   PetscCall(ISRestoreIndices(isrow, &irow));
345   PetscCall(ISRestoreIndices(iscol_local, &icol));
346   PetscCall(ISDestroy(&iscol_local));
347   *B = newmat;
348   PetscFunctionReturn(PETSC_SUCCESS);
349 }
350 
351 static PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A, PetscScalar **array)
352 {
353   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
354 
355   PetscFunctionBegin;
356   PetscCall(MatDenseRestoreArray(a->A, array));
357   PetscFunctionReturn(PETSC_SUCCESS);
358 }
359 
360 static PetscErrorCode MatDenseRestoreArrayRead_MPIDense(Mat A, PetscScalar **array)
361 {
362   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
363 
364   PetscFunctionBegin;
365   PetscCall(MatDenseRestoreArrayRead(a->A, (const PetscScalar **)array));
366   PetscFunctionReturn(PETSC_SUCCESS);
367 }
368 
369 static PetscErrorCode MatDenseRestoreArrayWrite_MPIDense(Mat A, PetscScalar **array)
370 {
371   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
372 
373   PetscFunctionBegin;
374   PetscCall(MatDenseRestoreArrayWrite(a->A, array));
375   PetscFunctionReturn(PETSC_SUCCESS);
376 }
377 
378 static PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat, MatAssemblyType mode)
379 {
380   Mat_MPIDense *mdn = (Mat_MPIDense *)mat->data;
381   PetscInt      nstash, reallocs;
382 
383   PetscFunctionBegin;
384   if (mdn->donotstash || mat->nooffprocentries) PetscFunctionReturn(PETSC_SUCCESS);
385 
386   PetscCall(MatStashScatterBegin_Private(mat, &mat->stash, mat->rmap->range));
387   PetscCall(MatStashGetInfo_Private(&mat->stash, &nstash, &reallocs));
388   PetscCall(PetscInfo(mdn->A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
389   PetscFunctionReturn(PETSC_SUCCESS);
390 }
391 
392 static PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat, MatAssemblyType mode)
393 {
394   Mat_MPIDense *mdn = (Mat_MPIDense *)mat->data;
395   PetscInt      i, *row, *col, flg, j, rstart, ncols;
396   PetscMPIInt   n;
397   PetscScalar  *val;
398 
399   PetscFunctionBegin;
400   if (!mdn->donotstash && !mat->nooffprocentries) {
401     /*  wait on receives */
402     while (1) {
403       PetscCall(MatStashScatterGetMesg_Private(&mat->stash, &n, &row, &col, &val, &flg));
404       if (!flg) break;
405 
406       for (i = 0; i < n;) {
407         /* Now identify the consecutive vals belonging to the same row */
408         for (j = i, rstart = row[j]; j < n; j++) {
409           if (row[j] != rstart) break;
410         }
411         if (j < n) ncols = j - i;
412         else ncols = n - i;
413         /* Now assemble all these values with a single function call */
414         PetscCall(MatSetValues_MPIDense(mat, 1, row + i, ncols, col + i, val + i, mat->insertmode));
415         i = j;
416       }
417     }
418     PetscCall(MatStashScatterEnd_Private(&mat->stash));
419   }
420 
421   PetscCall(MatAssemblyBegin(mdn->A, mode));
422   PetscCall(MatAssemblyEnd(mdn->A, mode));
423   PetscFunctionReturn(PETSC_SUCCESS);
424 }
425 
426 static PetscErrorCode MatZeroEntries_MPIDense(Mat A)
427 {
428   Mat_MPIDense *l = (Mat_MPIDense *)A->data;
429 
430   PetscFunctionBegin;
431   PetscCall(MatZeroEntries(l->A));
432   PetscFunctionReturn(PETSC_SUCCESS);
433 }
434 
435 static PetscErrorCode MatZeroRows_MPIDense(Mat A, PetscInt n, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
436 {
437   Mat_MPIDense *l = (Mat_MPIDense *)A->data;
438   PetscInt      i, len, *lrows;
439 
440   PetscFunctionBegin;
441   /* get locally owned rows */
442   PetscCall(PetscLayoutMapLocal(A->rmap, n, rows, &len, &lrows, NULL));
443   /* fix right-hand side if needed */
444   if (x && b) {
445     const PetscScalar *xx;
446     PetscScalar       *bb;
447 
448     PetscCall(VecGetArrayRead(x, &xx));
449     PetscCall(VecGetArrayWrite(b, &bb));
450     for (i = 0; i < len; ++i) bb[lrows[i]] = diag * xx[lrows[i]];
451     PetscCall(VecRestoreArrayRead(x, &xx));
452     PetscCall(VecRestoreArrayWrite(b, &bb));
453   }
454   PetscCall(MatZeroRows(l->A, len, lrows, 0.0, NULL, NULL));
455   if (diag != 0.0) {
456     Vec d;
457 
458     PetscCall(MatCreateVecs(A, NULL, &d));
459     PetscCall(VecSet(d, diag));
460     PetscCall(MatDiagonalSet(A, d, INSERT_VALUES));
461     PetscCall(VecDestroy(&d));
462   }
463   PetscCall(PetscFree(lrows));
464   PetscFunctionReturn(PETSC_SUCCESS);
465 }
466 
467 PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat, Vec, Vec);
468 PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat, Vec, Vec, Vec);
469 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat, Vec, Vec);
470 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat, Vec, Vec, Vec);
471 
472 static PetscErrorCode MatMult_MPIDense(Mat mat, Vec xx, Vec yy)
473 {
474   Mat_MPIDense      *mdn = (Mat_MPIDense *)mat->data;
475   const PetscScalar *ax;
476   PetscScalar       *ay;
477   PetscMemType       axmtype, aymtype;
478 
479   PetscFunctionBegin;
480   if (!mdn->Mvctx) PetscCall(MatSetUpMultiply_MPIDense(mat));
481   PetscCall(VecGetArrayReadAndMemType(xx, &ax, &axmtype));
482   PetscCall(VecGetArrayAndMemType(mdn->lvec, &ay, &aymtype));
483   PetscCall(PetscSFBcastWithMemTypeBegin(mdn->Mvctx, MPIU_SCALAR, axmtype, ax, aymtype, ay, MPI_REPLACE));
484   PetscCall(PetscSFBcastEnd(mdn->Mvctx, MPIU_SCALAR, ax, ay, MPI_REPLACE));
485   PetscCall(VecRestoreArrayAndMemType(mdn->lvec, &ay));
486   PetscCall(VecRestoreArrayReadAndMemType(xx, &ax));
487   PetscCall((*mdn->A->ops->mult)(mdn->A, mdn->lvec, yy));
488   PetscFunctionReturn(PETSC_SUCCESS);
489 }
490 
491 static PetscErrorCode MatMultAddColumnRange_MPIDense(Mat mat, Vec xx, Vec yy, Vec zz, PetscInt c_start, PetscInt c_end)
492 {
493   Mat_MPIDense      *mdn = (Mat_MPIDense *)mat->data;
494   const PetscScalar *ax;
495   PetscScalar       *ay;
496   PetscMemType       axmtype, aymtype;
497 
498   PetscFunctionBegin;
499   if (!mdn->Mvctx) PetscCall(MatSetUpMultiply_MPIDense(mat));
500   PetscCall(VecGetArrayReadAndMemType(xx, &ax, &axmtype));
501   PetscCall(VecGetArrayAndMemType(mdn->lvec, &ay, &aymtype));
502   PetscCall(PetscSFBcastWithMemTypeBegin(mdn->Mvctx, MPIU_SCALAR, axmtype, ax, aymtype, ay, MPI_REPLACE));
503   PetscCall(PetscSFBcastEnd(mdn->Mvctx, MPIU_SCALAR, ax, ay, MPI_REPLACE));
504   PetscCall(VecRestoreArrayAndMemType(mdn->lvec, &ay));
505   PetscCall(VecRestoreArrayReadAndMemType(xx, &ax));
506   PetscUseMethod(mdn->A, "MatMultAddColumnRange_C", (Mat, Vec, Vec, Vec, PetscInt, PetscInt), (mdn->A, mdn->lvec, yy, zz, c_start, c_end));
507   PetscFunctionReturn(PETSC_SUCCESS);
508 }
509 
510 static PetscErrorCode MatMultAdd_MPIDense(Mat mat, Vec xx, Vec yy, Vec zz)
511 {
512   Mat_MPIDense      *mdn = (Mat_MPIDense *)mat->data;
513   const PetscScalar *ax;
514   PetscScalar       *ay;
515   PetscMemType       axmtype, aymtype;
516 
517   PetscFunctionBegin;
518   if (!mdn->Mvctx) PetscCall(MatSetUpMultiply_MPIDense(mat));
519   PetscCall(VecGetArrayReadAndMemType(xx, &ax, &axmtype));
520   PetscCall(VecGetArrayAndMemType(mdn->lvec, &ay, &aymtype));
521   PetscCall(PetscSFBcastWithMemTypeBegin(mdn->Mvctx, MPIU_SCALAR, axmtype, ax, aymtype, ay, MPI_REPLACE));
522   PetscCall(PetscSFBcastEnd(mdn->Mvctx, MPIU_SCALAR, ax, ay, MPI_REPLACE));
523   PetscCall(VecRestoreArrayAndMemType(mdn->lvec, &ay));
524   PetscCall(VecRestoreArrayReadAndMemType(xx, &ax));
525   PetscCall((*mdn->A->ops->multadd)(mdn->A, mdn->lvec, yy, zz));
526   PetscFunctionReturn(PETSC_SUCCESS);
527 }
528 
529 static PetscErrorCode MatMultHermitianTransposeColumnRange_MPIDense(Mat A, Vec xx, Vec yy, PetscInt c_start, PetscInt c_end)
530 {
531   Mat_MPIDense      *a = (Mat_MPIDense *)A->data;
532   const PetscScalar *ax;
533   PetscScalar       *ay;
534   PetscMemType       axmtype, aymtype;
535 
536   PetscFunctionBegin;
537   if (!a->Mvctx) PetscCall(MatSetUpMultiply_MPIDense(A));
538   PetscCall(VecSet(yy, 0.0));
539   PetscUseMethod(a->A, "MatMultHermitianTransposeColumnRange_C", (Mat, Vec, Vec, PetscInt, PetscInt), (a->A, xx, a->lvec, c_start, c_end));
540   PetscCall(VecGetArrayReadAndMemType(a->lvec, &ax, &axmtype));
541   PetscCall(VecGetArrayAndMemType(yy, &ay, &aymtype));
542   PetscCall(PetscSFReduceWithMemTypeBegin(a->Mvctx, MPIU_SCALAR, axmtype, ax, aymtype, ay, MPIU_SUM));
543   PetscCall(PetscSFReduceEnd(a->Mvctx, MPIU_SCALAR, ax, ay, MPIU_SUM));
544   PetscCall(VecRestoreArrayReadAndMemType(a->lvec, &ax));
545   PetscCall(VecRestoreArrayAndMemType(yy, &ay));
546   PetscFunctionReturn(PETSC_SUCCESS);
547 }
548 
549 static PetscErrorCode MatMultTransposeKernel_MPIDense(Mat A, Vec xx, Vec yy, PetscBool herm)
550 {
551   Mat_MPIDense      *a = (Mat_MPIDense *)A->data;
552   const PetscScalar *ax;
553   PetscScalar       *ay;
554   PetscMemType       axmtype, aymtype;
555 
556   PetscFunctionBegin;
557   if (!a->Mvctx) PetscCall(MatSetUpMultiply_MPIDense(A));
558   PetscCall(VecSet(yy, 0.0));
559   if (herm) PetscCall((*a->A->ops->multhermitiantranspose)(a->A, xx, a->lvec));
560   else PetscCall((*a->A->ops->multtranspose)(a->A, xx, a->lvec));
561   PetscCall(VecGetArrayReadAndMemType(a->lvec, &ax, &axmtype));
562   PetscCall(VecGetArrayAndMemType(yy, &ay, &aymtype));
563   PetscCall(PetscSFReduceWithMemTypeBegin(a->Mvctx, MPIU_SCALAR, axmtype, ax, aymtype, ay, MPIU_SUM));
564   PetscCall(PetscSFReduceEnd(a->Mvctx, MPIU_SCALAR, ax, ay, MPIU_SUM));
565   PetscCall(VecRestoreArrayReadAndMemType(a->lvec, &ax));
566   PetscCall(VecRestoreArrayAndMemType(yy, &ay));
567   PetscFunctionReturn(PETSC_SUCCESS);
568 }
569 
570 static PetscErrorCode MatMultHermitianTransposeAddColumnRange_MPIDense(Mat A, Vec xx, Vec yy, Vec zz, PetscInt c_start, PetscInt c_end)
571 {
572   Mat_MPIDense      *a = (Mat_MPIDense *)A->data;
573   const PetscScalar *ax;
574   PetscScalar       *ay;
575   PetscMemType       axmtype, aymtype;
576 
577   PetscFunctionBegin;
578   if (!a->Mvctx) PetscCall(MatSetUpMultiply_MPIDense(A));
579   PetscCall(VecCopy(yy, zz));
580   PetscMPIInt rank;
581   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
582   PetscUseMethod(a->A, "MatMultHermitianTransposeColumnRange_C", (Mat, Vec, Vec, PetscInt, PetscInt), (a->A, xx, a->lvec, c_start, c_end));
583   PetscCall(VecGetArrayReadAndMemType(a->lvec, &ax, &axmtype));
584   PetscCall(VecGetArrayAndMemType(zz, &ay, &aymtype));
585   PetscCall(PetscSFReduceWithMemTypeBegin(a->Mvctx, MPIU_SCALAR, axmtype, ax, aymtype, ay, MPIU_SUM));
586   PetscCall(PetscSFReduceEnd(a->Mvctx, MPIU_SCALAR, ax, ay, MPIU_SUM));
587   PetscCall(VecRestoreArrayReadAndMemType(a->lvec, &ax));
588   PetscCall(VecRestoreArrayAndMemType(zz, &ay));
589   PetscFunctionReturn(PETSC_SUCCESS);
590 }
591 
592 static PetscErrorCode MatMultTransposeAddKernel_MPIDense(Mat A, Vec xx, Vec yy, Vec zz, PetscBool herm)
593 {
594   Mat_MPIDense      *a = (Mat_MPIDense *)A->data;
595   const PetscScalar *ax;
596   PetscScalar       *ay;
597   PetscMemType       axmtype, aymtype;
598 
599   PetscFunctionBegin;
600   if (!a->Mvctx) PetscCall(MatSetUpMultiply_MPIDense(A));
601   PetscCall(VecCopy(yy, zz));
602   if (herm) PetscCall((*a->A->ops->multhermitiantranspose)(a->A, xx, a->lvec));
603   else PetscCall((*a->A->ops->multtranspose)(a->A, xx, a->lvec));
604   PetscCall(VecGetArrayReadAndMemType(a->lvec, &ax, &axmtype));
605   PetscCall(VecGetArrayAndMemType(zz, &ay, &aymtype));
606   PetscCall(PetscSFReduceWithMemTypeBegin(a->Mvctx, MPIU_SCALAR, axmtype, ax, aymtype, ay, MPIU_SUM));
607   PetscCall(PetscSFReduceEnd(a->Mvctx, MPIU_SCALAR, ax, ay, MPIU_SUM));
608   PetscCall(VecRestoreArrayReadAndMemType(a->lvec, &ax));
609   PetscCall(VecRestoreArrayAndMemType(zz, &ay));
610   PetscFunctionReturn(PETSC_SUCCESS);
611 }
612 
613 static PetscErrorCode MatMultTranspose_MPIDense(Mat A, Vec xx, Vec yy)
614 {
615   PetscFunctionBegin;
616   PetscCall(MatMultTransposeKernel_MPIDense(A, xx, yy, PETSC_FALSE));
617   PetscFunctionReturn(PETSC_SUCCESS);
618 }
619 
620 static PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A, Vec xx, Vec yy, Vec zz)
621 {
622   PetscFunctionBegin;
623   PetscCall(MatMultTransposeAddKernel_MPIDense(A, xx, yy, zz, PETSC_FALSE));
624   PetscFunctionReturn(PETSC_SUCCESS);
625 }
626 
627 static PetscErrorCode MatMultHermitianTranspose_MPIDense(Mat A, Vec xx, Vec yy)
628 {
629   PetscFunctionBegin;
630   PetscCall(MatMultTransposeKernel_MPIDense(A, xx, yy, PETSC_TRUE));
631   PetscFunctionReturn(PETSC_SUCCESS);
632 }
633 
634 static PetscErrorCode MatMultHermitianTransposeAdd_MPIDense(Mat A, Vec xx, Vec yy, Vec zz)
635 {
636   PetscFunctionBegin;
637   PetscCall(MatMultTransposeAddKernel_MPIDense(A, xx, yy, zz, PETSC_TRUE));
638   PetscFunctionReturn(PETSC_SUCCESS);
639 }
640 
641 PetscErrorCode MatGetDiagonal_MPIDense(Mat A, Vec v)
642 {
643   Mat_MPIDense      *a = (Mat_MPIDense *)A->data;
644   PetscInt           lda, len, i, nl, ng, m = A->rmap->n, radd;
645   PetscScalar       *x;
646   const PetscScalar *av;
647 
648   PetscFunctionBegin;
649   PetscCall(VecGetArray(v, &x));
650   PetscCall(VecGetSize(v, &ng));
651   PetscCheck(ng == A->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming mat and vec");
652   PetscCall(VecGetLocalSize(v, &nl));
653   len  = PetscMin(a->A->rmap->n, a->A->cmap->n);
654   radd = A->rmap->rstart * m;
655   PetscCall(MatDenseGetArrayRead(a->A, &av));
656   PetscCall(MatDenseGetLDA(a->A, &lda));
657   for (i = 0; i < len; i++) x[i] = av[radd + i * lda + i];
658   PetscCall(MatDenseRestoreArrayRead(a->A, &av));
659   if (nl - i > 0) PetscCall(PetscArrayzero(x + i, nl - i));
660   PetscCall(VecRestoreArray(v, &x));
661   PetscFunctionReturn(PETSC_SUCCESS);
662 }
663 
664 static PetscErrorCode MatDestroy_MPIDense(Mat mat)
665 {
666   Mat_MPIDense *mdn = (Mat_MPIDense *)mat->data;
667 
668   PetscFunctionBegin;
669   PetscCall(PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT, mat->rmap->N, mat->cmap->N));
670   PetscCall(MatStashDestroy_Private(&mat->stash));
671   PetscCheck(!mdn->vecinuse, PetscObjectComm((PetscObject)mat), PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first");
672   PetscCheck(!mdn->matinuse, PetscObjectComm((PetscObject)mat), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
673   PetscCall(MatDestroy(&mdn->A));
674   PetscCall(VecDestroy(&mdn->lvec));
675   PetscCall(PetscSFDestroy(&mdn->Mvctx));
676   PetscCall(VecDestroy(&mdn->cvec));
677   PetscCall(MatDestroy(&mdn->cmat));
678 
679   PetscCall(PetscFree(mat->data));
680   PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL));
681 
682   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetLDA_C", NULL));
683   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseSetLDA_C", NULL));
684   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetArray_C", NULL));
685   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreArray_C", NULL));
686   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetArrayRead_C", NULL));
687   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreArrayRead_C", NULL));
688   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetArrayWrite_C", NULL));
689   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreArrayWrite_C", NULL));
690   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDensePlaceArray_C", NULL));
691   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseResetArray_C", NULL));
692   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseReplaceArray_C", NULL));
693   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpiaij_mpidense_C", NULL));
694   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_mpiaij_C", NULL));
695 #if defined(PETSC_HAVE_ELEMENTAL)
696   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_elemental_C", NULL));
697 #endif
698 #if defined(PETSC_HAVE_SCALAPACK)
699   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_scalapack_C", NULL));
700 #endif
701   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIDenseSetPreallocation_C", NULL));
702   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaij_mpidense_C", NULL));
703   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidense_mpiaij_C", NULL));
704 #if defined(PETSC_HAVE_CUDA)
705   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaijcusparse_mpidense_C", NULL));
706   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidense_mpiaijcusparse_C", NULL));
707   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_mpidensecuda_C", NULL));
708   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidensecuda_mpidense_C", NULL));
709   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaij_mpidensecuda_C", NULL));
710   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaijcusparse_mpidensecuda_C", NULL));
711   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidensecuda_mpiaij_C", NULL));
712   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidensecuda_mpiaijcusparse_C", NULL));
713   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDAGetArray_C", NULL));
714   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDAGetArrayRead_C", NULL));
715   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDAGetArrayWrite_C", NULL));
716   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDARestoreArray_C", NULL));
717   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDARestoreArrayRead_C", NULL));
718   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDARestoreArrayWrite_C", NULL));
719   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDAPlaceArray_C", NULL));
720   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDAResetArray_C", NULL));
721   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDAReplaceArray_C", NULL));
722   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseCUDASetPreallocation_C", NULL));
723 #endif
724 #if defined(PETSC_HAVE_HIP)
725   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaijhipsparse_mpidense_C", NULL));
726   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidense_mpiaijhipsparse_C", NULL));
727   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_mpidensehip_C", NULL));
728   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidensehip_mpidense_C", NULL));
729   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaij_mpidensehip_C", NULL));
730   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaijhipsparse_mpidensehip_C", NULL));
731   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidensehip_mpiaij_C", NULL));
732   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidensehip_mpiaijhipsparse_C", NULL));
733   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPGetArray_C", NULL));
734   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPGetArrayRead_C", NULL));
735   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPGetArrayWrite_C", NULL));
736   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPRestoreArray_C", NULL));
737   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPRestoreArrayRead_C", NULL));
738   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPRestoreArrayWrite_C", NULL));
739   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPPlaceArray_C", NULL));
740   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPResetArray_C", NULL));
741   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPReplaceArray_C", NULL));
742   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseHIPSetPreallocation_C", NULL));
743 #endif
744   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetColumn_C", NULL));
745   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreColumn_C", NULL));
746   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetColumnVec_C", NULL));
747   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreColumnVec_C", NULL));
748   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetColumnVecRead_C", NULL));
749   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreColumnVecRead_C", NULL));
750   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetColumnVecWrite_C", NULL));
751   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreColumnVecWrite_C", NULL));
752   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetSubMatrix_C", NULL));
753   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreSubMatrix_C", NULL));
754   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMultAddColumnRange_C", NULL));
755   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMultHermitianTransposeColumnRange_C", NULL));
756   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMultHermitianTransposeAddColumnRange_C", NULL));
757 
758   PetscCall(PetscObjectCompose((PetscObject)mat, "DiagonalBlock", NULL));
759   PetscFunctionReturn(PETSC_SUCCESS);
760 }
761 
762 #include <petscdraw.h>
763 static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat, PetscViewer viewer)
764 {
765   Mat_MPIDense     *mdn = (Mat_MPIDense *)mat->data;
766   PetscMPIInt       rank;
767   PetscViewerType   vtype;
768   PetscBool         iascii, isdraw;
769   PetscViewer       sviewer;
770   PetscViewerFormat format;
771 
772   PetscFunctionBegin;
773   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat), &rank));
774   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
775   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
776   if (iascii) {
777     PetscCall(PetscViewerGetType(viewer, &vtype));
778     PetscCall(PetscViewerGetFormat(viewer, &format));
779     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
780       MatInfo info;
781       PetscCall(MatGetInfo(mat, MAT_LOCAL, &info));
782       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
783       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] local rows %" PetscInt_FMT " nz %" PetscInt_FMT " nz alloced %" PetscInt_FMT " mem %" PetscInt_FMT " \n", rank, mat->rmap->n, (PetscInt)info.nz_used, (PetscInt)info.nz_allocated,
784                                                    (PetscInt)info.memory));
785       PetscCall(PetscViewerFlush(viewer));
786       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
787       if (mdn->Mvctx) PetscCall(PetscSFView(mdn->Mvctx, viewer));
788       PetscFunctionReturn(PETSC_SUCCESS);
789     } else if (format == PETSC_VIEWER_ASCII_INFO) {
790       PetscFunctionReturn(PETSC_SUCCESS);
791     }
792   } else if (isdraw) {
793     PetscDraw draw;
794     PetscBool isnull;
795 
796     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
797     PetscCall(PetscDrawIsNull(draw, &isnull));
798     if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
799   }
800 
801   {
802     /* assemble the entire matrix onto first processor. */
803     Mat          A;
804     PetscInt     M = mat->rmap->N, N = mat->cmap->N, m, row, i, nz;
805     PetscInt    *cols;
806     PetscScalar *vals;
807 
808     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &A));
809     if (rank == 0) {
810       PetscCall(MatSetSizes(A, M, N, M, N));
811     } else {
812       PetscCall(MatSetSizes(A, 0, 0, M, N));
813     }
814     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
815     PetscCall(MatSetType(A, MATMPIDENSE));
816     PetscCall(MatMPIDenseSetPreallocation(A, NULL));
817 
818     /* Copy the matrix ... This isn't the most efficient means,
819        but it's quick for now */
820     A->insertmode = INSERT_VALUES;
821 
822     row = mat->rmap->rstart;
823     m   = mdn->A->rmap->n;
824     for (i = 0; i < m; i++) {
825       PetscCall(MatGetRow_MPIDense(mat, row, &nz, &cols, &vals));
826       PetscCall(MatSetValues_MPIDense(A, 1, &row, nz, cols, vals, INSERT_VALUES));
827       PetscCall(MatRestoreRow_MPIDense(mat, row, &nz, &cols, &vals));
828       row++;
829     }
830 
831     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
832     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
833     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
834     if (rank == 0) {
835       PetscCall(PetscObjectSetName((PetscObject)((Mat_MPIDense *)A->data)->A, ((PetscObject)mat)->name));
836       PetscCall(MatView_SeqDense(((Mat_MPIDense *)A->data)->A, sviewer));
837     }
838     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
839     PetscCall(MatDestroy(&A));
840   }
841   PetscFunctionReturn(PETSC_SUCCESS);
842 }
843 
844 static PetscErrorCode MatView_MPIDense(Mat mat, PetscViewer viewer)
845 {
846   PetscBool iascii, isbinary, isdraw, issocket;
847 
848   PetscFunctionBegin;
849   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
850   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
851   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSOCKET, &issocket));
852   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
853 
854   if (iascii || issocket || isdraw) {
855     PetscCall(MatView_MPIDense_ASCIIorDraworSocket(mat, viewer));
856   } else if (isbinary) PetscCall(MatView_Dense_Binary(mat, viewer));
857   PetscFunctionReturn(PETSC_SUCCESS);
858 }
859 
860 static PetscErrorCode MatGetInfo_MPIDense(Mat A, MatInfoType flag, MatInfo *info)
861 {
862   Mat_MPIDense  *mat = (Mat_MPIDense *)A->data;
863   Mat            mdn = mat->A;
864   PetscLogDouble isend[5], irecv[5];
865 
866   PetscFunctionBegin;
867   info->block_size = 1.0;
868 
869   PetscCall(MatGetInfo(mdn, MAT_LOCAL, info));
870 
871   isend[0] = info->nz_used;
872   isend[1] = info->nz_allocated;
873   isend[2] = info->nz_unneeded;
874   isend[3] = info->memory;
875   isend[4] = info->mallocs;
876   if (flag == MAT_LOCAL) {
877     info->nz_used      = isend[0];
878     info->nz_allocated = isend[1];
879     info->nz_unneeded  = isend[2];
880     info->memory       = isend[3];
881     info->mallocs      = isend[4];
882   } else if (flag == MAT_GLOBAL_MAX) {
883     PetscCall(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_MAX, PetscObjectComm((PetscObject)A)));
884 
885     info->nz_used      = irecv[0];
886     info->nz_allocated = irecv[1];
887     info->nz_unneeded  = irecv[2];
888     info->memory       = irecv[3];
889     info->mallocs      = irecv[4];
890   } else if (flag == MAT_GLOBAL_SUM) {
891     PetscCall(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_SUM, PetscObjectComm((PetscObject)A)));
892 
893     info->nz_used      = irecv[0];
894     info->nz_allocated = irecv[1];
895     info->nz_unneeded  = irecv[2];
896     info->memory       = irecv[3];
897     info->mallocs      = irecv[4];
898   }
899   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
900   info->fill_ratio_needed = 0;
901   info->factor_mallocs    = 0;
902   PetscFunctionReturn(PETSC_SUCCESS);
903 }
904 
905 static PetscErrorCode MatSetOption_MPIDense(Mat A, MatOption op, PetscBool flg)
906 {
907   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
908 
909   PetscFunctionBegin;
910   switch (op) {
911   case MAT_NEW_NONZERO_LOCATIONS:
912   case MAT_NEW_NONZERO_LOCATION_ERR:
913   case MAT_NEW_NONZERO_ALLOCATION_ERR:
914     MatCheckPreallocated(A, 1);
915     PetscCall(MatSetOption(a->A, op, flg));
916     break;
917   case MAT_ROW_ORIENTED:
918     MatCheckPreallocated(A, 1);
919     a->roworiented = flg;
920     PetscCall(MatSetOption(a->A, op, flg));
921     break;
922   case MAT_FORCE_DIAGONAL_ENTRIES:
923   case MAT_KEEP_NONZERO_PATTERN:
924   case MAT_USE_HASH_TABLE:
925   case MAT_SORTED_FULL:
926     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
927     break;
928   case MAT_IGNORE_OFF_PROC_ENTRIES:
929     a->donotstash = flg;
930     break;
931   case MAT_SYMMETRIC:
932   case MAT_STRUCTURALLY_SYMMETRIC:
933   case MAT_HERMITIAN:
934   case MAT_SYMMETRY_ETERNAL:
935   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
936   case MAT_SPD:
937   case MAT_IGNORE_LOWER_TRIANGULAR:
938   case MAT_IGNORE_ZERO_ENTRIES:
939   case MAT_SPD_ETERNAL:
940     /* if the diagonal matrix is square it inherits some of the properties above */
941     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
942     break;
943   default:
944     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %s", MatOptions[op]);
945   }
946   PetscFunctionReturn(PETSC_SUCCESS);
947 }
948 
949 static PetscErrorCode MatDiagonalScale_MPIDense(Mat A, Vec ll, Vec rr)
950 {
951   Mat_MPIDense      *mdn = (Mat_MPIDense *)A->data;
952   const PetscScalar *l;
953   PetscScalar        x, *v, *vv, *r;
954   PetscInt           i, j, s2a, s3a, s2, s3, m = mdn->A->rmap->n, n = mdn->A->cmap->n, lda;
955 
956   PetscFunctionBegin;
957   PetscCall(MatDenseGetArray(mdn->A, &vv));
958   PetscCall(MatDenseGetLDA(mdn->A, &lda));
959   PetscCall(MatGetLocalSize(A, &s2, &s3));
960   if (ll) {
961     PetscCall(VecGetLocalSize(ll, &s2a));
962     PetscCheck(s2a == s2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Left scaling vector non-conforming local size, %" PetscInt_FMT " != %" PetscInt_FMT, s2a, s2);
963     PetscCall(VecGetArrayRead(ll, &l));
964     for (i = 0; i < m; i++) {
965       x = l[i];
966       v = vv + i;
967       for (j = 0; j < n; j++) {
968         (*v) *= x;
969         v += lda;
970       }
971     }
972     PetscCall(VecRestoreArrayRead(ll, &l));
973     PetscCall(PetscLogFlops(1.0 * n * m));
974   }
975   if (rr) {
976     const PetscScalar *ar;
977 
978     PetscCall(VecGetLocalSize(rr, &s3a));
979     PetscCheck(s3a == s3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Right scaling vec non-conforming local size, %" PetscInt_FMT " != %" PetscInt_FMT ".", s3a, s3);
980     PetscCall(VecGetArrayRead(rr, &ar));
981     if (!mdn->Mvctx) PetscCall(MatSetUpMultiply_MPIDense(A));
982     PetscCall(VecGetArray(mdn->lvec, &r));
983     PetscCall(PetscSFBcastBegin(mdn->Mvctx, MPIU_SCALAR, ar, r, MPI_REPLACE));
984     PetscCall(PetscSFBcastEnd(mdn->Mvctx, MPIU_SCALAR, ar, r, MPI_REPLACE));
985     PetscCall(VecRestoreArrayRead(rr, &ar));
986     for (i = 0; i < n; i++) {
987       x = r[i];
988       v = vv + i * lda;
989       for (j = 0; j < m; j++) (*v++) *= x;
990     }
991     PetscCall(VecRestoreArray(mdn->lvec, &r));
992     PetscCall(PetscLogFlops(1.0 * n * m));
993   }
994   PetscCall(MatDenseRestoreArray(mdn->A, &vv));
995   PetscFunctionReturn(PETSC_SUCCESS);
996 }
997 
998 static PetscErrorCode MatNorm_MPIDense(Mat A, NormType type, PetscReal *nrm)
999 {
1000   Mat_MPIDense      *mdn = (Mat_MPIDense *)A->data;
1001   PetscInt           i, j;
1002   PetscMPIInt        size;
1003   PetscReal          sum = 0.0;
1004   const PetscScalar *av, *v;
1005 
1006   PetscFunctionBegin;
1007   PetscCall(MatDenseGetArrayRead(mdn->A, &av));
1008   v = av;
1009   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1010   if (size == 1) {
1011     PetscCall(MatNorm(mdn->A, type, nrm));
1012   } else {
1013     if (type == NORM_FROBENIUS) {
1014       for (i = 0; i < mdn->A->cmap->n * mdn->A->rmap->n; i++) {
1015         sum += PetscRealPart(PetscConj(*v) * (*v));
1016         v++;
1017       }
1018       PetscCall(MPIU_Allreduce(&sum, nrm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)A)));
1019       *nrm = PetscSqrtReal(*nrm);
1020       PetscCall(PetscLogFlops(2.0 * mdn->A->cmap->n * mdn->A->rmap->n));
1021     } else if (type == NORM_1) {
1022       PetscReal *tmp, *tmp2;
1023       PetscCall(PetscCalloc2(A->cmap->N, &tmp, A->cmap->N, &tmp2));
1024       *nrm = 0.0;
1025       v    = av;
1026       for (j = 0; j < mdn->A->cmap->n; j++) {
1027         for (i = 0; i < mdn->A->rmap->n; i++) {
1028           tmp[j] += PetscAbsScalar(*v);
1029           v++;
1030         }
1031       }
1032       PetscCall(MPIU_Allreduce(tmp, tmp2, A->cmap->N, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)A)));
1033       for (j = 0; j < A->cmap->N; j++) {
1034         if (tmp2[j] > *nrm) *nrm = tmp2[j];
1035       }
1036       PetscCall(PetscFree2(tmp, tmp2));
1037       PetscCall(PetscLogFlops(A->cmap->n * A->rmap->n));
1038     } else if (type == NORM_INFINITY) { /* max row norm */
1039       PetscReal ntemp;
1040       PetscCall(MatNorm(mdn->A, type, &ntemp));
1041       PetscCall(MPIU_Allreduce(&ntemp, nrm, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)A)));
1042     } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for two norm");
1043   }
1044   PetscCall(MatDenseRestoreArrayRead(mdn->A, &av));
1045   PetscFunctionReturn(PETSC_SUCCESS);
1046 }
1047 
1048 static PetscErrorCode MatTranspose_MPIDense(Mat A, MatReuse reuse, Mat *matout)
1049 {
1050   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1051   Mat           B;
1052   PetscInt      M = A->rmap->N, N = A->cmap->N, m, n, *rwork, rstart = A->rmap->rstart;
1053   PetscInt      j, i, lda;
1054   PetscScalar  *v;
1055 
1056   PetscFunctionBegin;
1057   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *matout));
1058   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
1059     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
1060     PetscCall(MatSetSizes(B, A->cmap->n, A->rmap->n, N, M));
1061     PetscCall(MatSetType(B, ((PetscObject)A)->type_name));
1062     PetscCall(MatMPIDenseSetPreallocation(B, NULL));
1063   } else B = *matout;
1064 
1065   m = a->A->rmap->n;
1066   n = a->A->cmap->n;
1067   PetscCall(MatDenseGetArrayRead(a->A, (const PetscScalar **)&v));
1068   PetscCall(MatDenseGetLDA(a->A, &lda));
1069   PetscCall(PetscMalloc1(m, &rwork));
1070   for (i = 0; i < m; i++) rwork[i] = rstart + i;
1071   for (j = 0; j < n; j++) {
1072     PetscCall(MatSetValues(B, 1, &j, m, rwork, v, INSERT_VALUES));
1073     v = PetscSafePointerPlusOffset(v, lda);
1074   }
1075   PetscCall(MatDenseRestoreArrayRead(a->A, (const PetscScalar **)&v));
1076   PetscCall(PetscFree(rwork));
1077   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1078   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1079   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
1080     *matout = B;
1081   } else {
1082     PetscCall(MatHeaderMerge(A, &B));
1083   }
1084   PetscFunctionReturn(PETSC_SUCCESS);
1085 }
1086 
1087 static PetscErrorCode       MatDuplicate_MPIDense(Mat, MatDuplicateOption, Mat *);
1088 PETSC_INTERN PetscErrorCode MatScale_MPIDense(Mat, PetscScalar);
1089 
1090 static PetscErrorCode MatSetUp_MPIDense(Mat A)
1091 {
1092   PetscFunctionBegin;
1093   PetscCall(PetscLayoutSetUp(A->rmap));
1094   PetscCall(PetscLayoutSetUp(A->cmap));
1095   if (!A->preallocated) PetscCall(MatMPIDenseSetPreallocation(A, NULL));
1096   PetscFunctionReturn(PETSC_SUCCESS);
1097 }
1098 
1099 static PetscErrorCode MatAXPY_MPIDense(Mat Y, PetscScalar alpha, Mat X, MatStructure str)
1100 {
1101   Mat_MPIDense *A = (Mat_MPIDense *)Y->data, *B = (Mat_MPIDense *)X->data;
1102 
1103   PetscFunctionBegin;
1104   PetscCall(MatAXPY(A->A, alpha, B->A, str));
1105   PetscFunctionReturn(PETSC_SUCCESS);
1106 }
1107 
1108 static PetscErrorCode MatConjugate_MPIDense(Mat mat)
1109 {
1110   Mat_MPIDense *a = (Mat_MPIDense *)mat->data;
1111 
1112   PetscFunctionBegin;
1113   PetscCall(MatConjugate(a->A));
1114   PetscFunctionReturn(PETSC_SUCCESS);
1115 }
1116 
1117 static PetscErrorCode MatRealPart_MPIDense(Mat A)
1118 {
1119   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1120 
1121   PetscFunctionBegin;
1122   PetscCall(MatRealPart(a->A));
1123   PetscFunctionReturn(PETSC_SUCCESS);
1124 }
1125 
1126 static PetscErrorCode MatImaginaryPart_MPIDense(Mat A)
1127 {
1128   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1129 
1130   PetscFunctionBegin;
1131   PetscCall(MatImaginaryPart(a->A));
1132   PetscFunctionReturn(PETSC_SUCCESS);
1133 }
1134 
1135 static PetscErrorCode MatGetColumnVector_MPIDense(Mat A, Vec v, PetscInt col)
1136 {
1137   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1138 
1139   PetscFunctionBegin;
1140   PetscCheck(a->A, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Missing local matrix");
1141   PetscCheck(a->A->ops->getcolumnvector, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Missing get column operation");
1142   PetscCall((*a->A->ops->getcolumnvector)(a->A, v, col));
1143   PetscFunctionReturn(PETSC_SUCCESS);
1144 }
1145 
1146 PETSC_INTERN PetscErrorCode MatGetColumnReductions_SeqDense(Mat, PetscInt, PetscReal *);
1147 
1148 static PetscErrorCode MatGetColumnReductions_MPIDense(Mat A, PetscInt type, PetscReal *reductions)
1149 {
1150   PetscInt      i, m, n;
1151   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1152   PetscReal    *work;
1153 
1154   PetscFunctionBegin;
1155   PetscCall(MatGetSize(A, &m, &n));
1156   PetscCall(PetscMalloc1(n, &work));
1157   if (type == REDUCTION_MEAN_REALPART) {
1158     PetscCall(MatGetColumnReductions_SeqDense(a->A, (PetscInt)REDUCTION_SUM_REALPART, work));
1159   } else if (type == REDUCTION_MEAN_IMAGINARYPART) {
1160     PetscCall(MatGetColumnReductions_SeqDense(a->A, (PetscInt)REDUCTION_SUM_IMAGINARYPART, work));
1161   } else {
1162     PetscCall(MatGetColumnReductions_SeqDense(a->A, type, work));
1163   }
1164   if (type == NORM_2) {
1165     for (i = 0; i < n; i++) work[i] *= work[i];
1166   }
1167   if (type == NORM_INFINITY) {
1168     PetscCall(MPIU_Allreduce(work, reductions, n, MPIU_REAL, MPIU_MAX, A->hdr.comm));
1169   } else {
1170     PetscCall(MPIU_Allreduce(work, reductions, n, MPIU_REAL, MPIU_SUM, A->hdr.comm));
1171   }
1172   PetscCall(PetscFree(work));
1173   if (type == NORM_2) {
1174     for (i = 0; i < n; i++) reductions[i] = PetscSqrtReal(reductions[i]);
1175   } else if (type == REDUCTION_MEAN_REALPART || type == REDUCTION_MEAN_IMAGINARYPART) {
1176     for (i = 0; i < n; i++) reductions[i] /= m;
1177   }
1178   PetscFunctionReturn(PETSC_SUCCESS);
1179 }
1180 
1181 static PetscErrorCode MatSetRandom_MPIDense(Mat x, PetscRandom rctx)
1182 {
1183   Mat_MPIDense *d = (Mat_MPIDense *)x->data;
1184 
1185   PetscFunctionBegin;
1186   PetscCall(MatSetRandom(d->A, rctx));
1187 #if defined(PETSC_HAVE_DEVICE)
1188   x->offloadmask = d->A->offloadmask;
1189 #endif
1190   PetscFunctionReturn(PETSC_SUCCESS);
1191 }
1192 
1193 static PetscErrorCode MatMissingDiagonal_MPIDense(Mat A, PetscBool *missing, PetscInt *d)
1194 {
1195   PetscFunctionBegin;
1196   *missing = PETSC_FALSE;
1197   PetscFunctionReturn(PETSC_SUCCESS);
1198 }
1199 
1200 static PetscErrorCode MatMatTransposeMultSymbolic_MPIDense_MPIDense(Mat, Mat, PetscReal, Mat);
1201 static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense(Mat, Mat, Mat);
1202 static PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat, Mat, PetscReal, Mat);
1203 static PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat, Mat, Mat);
1204 static PetscErrorCode MatEqual_MPIDense(Mat, Mat, PetscBool *);
1205 static PetscErrorCode MatLoad_MPIDense(Mat, PetscViewer);
1206 static PetscErrorCode MatProductSetFromOptions_MPIDense(Mat);
1207 
1208 static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
1209                                        MatGetRow_MPIDense,
1210                                        MatRestoreRow_MPIDense,
1211                                        MatMult_MPIDense,
1212                                        /*  4*/ MatMultAdd_MPIDense,
1213                                        MatMultTranspose_MPIDense,
1214                                        MatMultTransposeAdd_MPIDense,
1215                                        NULL,
1216                                        NULL,
1217                                        NULL,
1218                                        /* 10*/ NULL,
1219                                        NULL,
1220                                        NULL,
1221                                        NULL,
1222                                        MatTranspose_MPIDense,
1223                                        /* 15*/ MatGetInfo_MPIDense,
1224                                        MatEqual_MPIDense,
1225                                        MatGetDiagonal_MPIDense,
1226                                        MatDiagonalScale_MPIDense,
1227                                        MatNorm_MPIDense,
1228                                        /* 20*/ MatAssemblyBegin_MPIDense,
1229                                        MatAssemblyEnd_MPIDense,
1230                                        MatSetOption_MPIDense,
1231                                        MatZeroEntries_MPIDense,
1232                                        /* 24*/ MatZeroRows_MPIDense,
1233                                        NULL,
1234                                        NULL,
1235                                        NULL,
1236                                        NULL,
1237                                        /* 29*/ MatSetUp_MPIDense,
1238                                        NULL,
1239                                        NULL,
1240                                        MatGetDiagonalBlock_MPIDense,
1241                                        NULL,
1242                                        /* 34*/ MatDuplicate_MPIDense,
1243                                        NULL,
1244                                        NULL,
1245                                        NULL,
1246                                        NULL,
1247                                        /* 39*/ MatAXPY_MPIDense,
1248                                        MatCreateSubMatrices_MPIDense,
1249                                        NULL,
1250                                        MatGetValues_MPIDense,
1251                                        MatCopy_MPIDense,
1252                                        /* 44*/ NULL,
1253                                        MatScale_MPIDense,
1254                                        MatShift_MPIDense,
1255                                        NULL,
1256                                        NULL,
1257                                        /* 49*/ MatSetRandom_MPIDense,
1258                                        NULL,
1259                                        NULL,
1260                                        NULL,
1261                                        NULL,
1262                                        /* 54*/ NULL,
1263                                        NULL,
1264                                        NULL,
1265                                        NULL,
1266                                        NULL,
1267                                        /* 59*/ MatCreateSubMatrix_MPIDense,
1268                                        MatDestroy_MPIDense,
1269                                        MatView_MPIDense,
1270                                        NULL,
1271                                        NULL,
1272                                        /* 64*/ NULL,
1273                                        NULL,
1274                                        NULL,
1275                                        NULL,
1276                                        NULL,
1277                                        /* 69*/ NULL,
1278                                        NULL,
1279                                        NULL,
1280                                        NULL,
1281                                        NULL,
1282                                        /* 74*/ NULL,
1283                                        NULL,
1284                                        NULL,
1285                                        NULL,
1286                                        NULL,
1287                                        /* 79*/ NULL,
1288                                        NULL,
1289                                        NULL,
1290                                        NULL,
1291                                        /* 83*/ MatLoad_MPIDense,
1292                                        NULL,
1293                                        NULL,
1294                                        NULL,
1295                                        NULL,
1296                                        NULL,
1297                                        /* 89*/ NULL,
1298                                        NULL,
1299                                        NULL,
1300                                        NULL,
1301                                        NULL,
1302                                        /* 94*/ NULL,
1303                                        NULL,
1304                                        MatMatTransposeMultSymbolic_MPIDense_MPIDense,
1305                                        MatMatTransposeMultNumeric_MPIDense_MPIDense,
1306                                        NULL,
1307                                        /* 99*/ MatProductSetFromOptions_MPIDense,
1308                                        NULL,
1309                                        NULL,
1310                                        MatConjugate_MPIDense,
1311                                        NULL,
1312                                        /*104*/ NULL,
1313                                        MatRealPart_MPIDense,
1314                                        MatImaginaryPart_MPIDense,
1315                                        NULL,
1316                                        NULL,
1317                                        /*109*/ NULL,
1318                                        NULL,
1319                                        NULL,
1320                                        MatGetColumnVector_MPIDense,
1321                                        MatMissingDiagonal_MPIDense,
1322                                        /*114*/ NULL,
1323                                        NULL,
1324                                        NULL,
1325                                        NULL,
1326                                        NULL,
1327                                        /*119*/ NULL,
1328                                        NULL,
1329                                        MatMultHermitianTranspose_MPIDense,
1330                                        MatMultHermitianTransposeAdd_MPIDense,
1331                                        NULL,
1332                                        /*124*/ NULL,
1333                                        MatGetColumnReductions_MPIDense,
1334                                        NULL,
1335                                        NULL,
1336                                        NULL,
1337                                        /*129*/ NULL,
1338                                        NULL,
1339                                        MatTransposeMatMultSymbolic_MPIDense_MPIDense,
1340                                        MatTransposeMatMultNumeric_MPIDense_MPIDense,
1341                                        NULL,
1342                                        /*134*/ NULL,
1343                                        NULL,
1344                                        NULL,
1345                                        NULL,
1346                                        NULL,
1347                                        /*139*/ NULL,
1348                                        NULL,
1349                                        NULL,
1350                                        NULL,
1351                                        NULL,
1352                                        MatCreateMPIMatConcatenateSeqMat_MPIDense,
1353                                        /*145*/ NULL,
1354                                        NULL,
1355                                        NULL,
1356                                        NULL,
1357                                        NULL,
1358                                        /*150*/ NULL,
1359                                        NULL,
1360                                        NULL,
1361                                        NULL,
1362                                        NULL,
1363                                        NULL};
1364 
1365 static PetscErrorCode MatMPIDenseSetPreallocation_MPIDense(Mat mat, PetscScalar *data)
1366 {
1367   Mat_MPIDense *a     = (Mat_MPIDense *)mat->data;
1368   MatType       mtype = MATSEQDENSE;
1369 
1370   PetscFunctionBegin;
1371   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)mat), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
1372   PetscCall(PetscLayoutSetUp(mat->rmap));
1373   PetscCall(PetscLayoutSetUp(mat->cmap));
1374   if (!a->A) {
1375     PetscCall(MatCreate(PETSC_COMM_SELF, &a->A));
1376     PetscCall(MatSetSizes(a->A, mat->rmap->n, mat->cmap->N, mat->rmap->n, mat->cmap->N));
1377   }
1378 #if defined(PETSC_HAVE_CUDA)
1379   PetscBool iscuda;
1380   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIDENSECUDA, &iscuda));
1381   if (iscuda) mtype = MATSEQDENSECUDA;
1382 #endif
1383 #if defined(PETSC_HAVE_HIP)
1384   PetscBool iship;
1385   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIDENSEHIP, &iship));
1386   if (iship) mtype = MATSEQDENSEHIP;
1387 #endif
1388   PetscCall(MatSetType(a->A, mtype));
1389   PetscCall(MatSeqDenseSetPreallocation(a->A, data));
1390 #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
1391   mat->offloadmask = a->A->offloadmask;
1392 #endif
1393   mat->preallocated = PETSC_TRUE;
1394   mat->assembled    = PETSC_TRUE;
1395   PetscFunctionReturn(PETSC_SUCCESS);
1396 }
1397 
1398 PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPIDense(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1399 {
1400   Mat B, C;
1401 
1402   PetscFunctionBegin;
1403   PetscCall(MatMPIAIJGetLocalMat(A, MAT_INITIAL_MATRIX, &C));
1404   PetscCall(MatConvert_SeqAIJ_SeqDense(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &B));
1405   PetscCall(MatDestroy(&C));
1406   if (reuse == MAT_REUSE_MATRIX) {
1407     C = *newmat;
1408   } else C = NULL;
1409   PetscCall(MatCreateMPIMatConcatenateSeqMat(PetscObjectComm((PetscObject)A), B, A->cmap->n, !C ? MAT_INITIAL_MATRIX : MAT_REUSE_MATRIX, &C));
1410   PetscCall(MatDestroy(&B));
1411   if (reuse == MAT_INPLACE_MATRIX) {
1412     PetscCall(MatHeaderReplace(A, &C));
1413   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = C;
1414   PetscFunctionReturn(PETSC_SUCCESS);
1415 }
1416 
1417 static PetscErrorCode MatConvert_MPIDense_MPIAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1418 {
1419   Mat B, C;
1420 
1421   PetscFunctionBegin;
1422   PetscCall(MatDenseGetLocalMatrix(A, &C));
1423   PetscCall(MatConvert_SeqDense_SeqAIJ(C, MATSEQAIJ, MAT_INITIAL_MATRIX, &B));
1424   if (reuse == MAT_REUSE_MATRIX) {
1425     C = *newmat;
1426   } else C = NULL;
1427   PetscCall(MatCreateMPIMatConcatenateSeqMat(PetscObjectComm((PetscObject)A), B, A->cmap->n, !C ? MAT_INITIAL_MATRIX : MAT_REUSE_MATRIX, &C));
1428   PetscCall(MatDestroy(&B));
1429   if (reuse == MAT_INPLACE_MATRIX) {
1430     PetscCall(MatHeaderReplace(A, &C));
1431   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = C;
1432   PetscFunctionReturn(PETSC_SUCCESS);
1433 }
1434 
1435 #if defined(PETSC_HAVE_ELEMENTAL)
1436 PETSC_INTERN PetscErrorCode MatConvert_MPIDense_Elemental(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1437 {
1438   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1439   Mat           mat_elemental;
1440   PetscScalar  *v;
1441   PetscInt      m = A->rmap->n, N = A->cmap->N, rstart = A->rmap->rstart, i, *rows, *cols, lda;
1442 
1443   PetscFunctionBegin;
1444   if (reuse == MAT_REUSE_MATRIX) {
1445     mat_elemental = *newmat;
1446     PetscCall(MatZeroEntries(*newmat));
1447   } else {
1448     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental));
1449     PetscCall(MatSetSizes(mat_elemental, PETSC_DECIDE, PETSC_DECIDE, A->rmap->N, A->cmap->N));
1450     PetscCall(MatSetType(mat_elemental, MATELEMENTAL));
1451     PetscCall(MatSetUp(mat_elemental));
1452     PetscCall(MatSetOption(mat_elemental, MAT_ROW_ORIENTED, PETSC_FALSE));
1453   }
1454 
1455   PetscCall(PetscMalloc2(m, &rows, N, &cols));
1456   for (i = 0; i < N; i++) cols[i] = i;
1457   for (i = 0; i < m; i++) rows[i] = rstart + i;
1458 
1459   /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
1460   PetscCall(MatDenseGetArray(A, &v));
1461   PetscCall(MatDenseGetLDA(a->A, &lda));
1462   if (lda == m) PetscCall(MatSetValues(mat_elemental, m, rows, N, cols, v, ADD_VALUES));
1463   else {
1464     for (i = 0; i < N; i++) PetscCall(MatSetValues(mat_elemental, m, rows, 1, &i, v + lda * i, ADD_VALUES));
1465   }
1466   PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY));
1467   PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY));
1468   PetscCall(MatDenseRestoreArray(A, &v));
1469   PetscCall(PetscFree2(rows, cols));
1470 
1471   if (reuse == MAT_INPLACE_MATRIX) {
1472     PetscCall(MatHeaderReplace(A, &mat_elemental));
1473   } else {
1474     *newmat = mat_elemental;
1475   }
1476   PetscFunctionReturn(PETSC_SUCCESS);
1477 }
1478 #endif
1479 
1480 static PetscErrorCode MatDenseGetColumn_MPIDense(Mat A, PetscInt col, PetscScalar **vals)
1481 {
1482   Mat_MPIDense *mat = (Mat_MPIDense *)A->data;
1483 
1484   PetscFunctionBegin;
1485   PetscCall(MatDenseGetColumn(mat->A, col, vals));
1486   PetscFunctionReturn(PETSC_SUCCESS);
1487 }
1488 
1489 static PetscErrorCode MatDenseRestoreColumn_MPIDense(Mat A, PetscScalar **vals)
1490 {
1491   Mat_MPIDense *mat = (Mat_MPIDense *)A->data;
1492 
1493   PetscFunctionBegin;
1494   PetscCall(MatDenseRestoreColumn(mat->A, vals));
1495   PetscFunctionReturn(PETSC_SUCCESS);
1496 }
1497 
1498 PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm comm, Mat inmat, PetscInt n, MatReuse scall, Mat *outmat)
1499 {
1500   Mat_MPIDense *mat;
1501   PetscInt      m, nloc, N;
1502 
1503   PetscFunctionBegin;
1504   PetscCall(MatGetSize(inmat, &m, &N));
1505   PetscCall(MatGetLocalSize(inmat, NULL, &nloc));
1506   if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
1507     PetscInt sum;
1508 
1509     if (n == PETSC_DECIDE) PetscCall(PetscSplitOwnership(comm, &n, &N));
1510     /* Check sum(n) = N */
1511     PetscCall(MPIU_Allreduce(&n, &sum, 1, MPIU_INT, MPI_SUM, comm));
1512     PetscCheck(sum == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Sum of local columns %" PetscInt_FMT " != global columns %" PetscInt_FMT, sum, N);
1513 
1514     PetscCall(MatCreateDense(comm, m, n, PETSC_DETERMINE, N, NULL, outmat));
1515     PetscCall(MatSetOption(*outmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1516   }
1517 
1518   /* numeric phase */
1519   mat = (Mat_MPIDense *)(*outmat)->data;
1520   PetscCall(MatCopy(inmat, mat->A, SAME_NONZERO_PATTERN));
1521   PetscFunctionReturn(PETSC_SUCCESS);
1522 }
1523 
1524 PetscErrorCode MatDenseGetColumnVec_MPIDense(Mat A, PetscInt col, Vec *v)
1525 {
1526   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1527   PetscInt      lda;
1528 
1529   PetscFunctionBegin;
1530   PetscCheck(!a->vecinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first");
1531   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
1532   if (!a->cvec) PetscCall(MatDenseCreateColumnVec_Private(A, &a->cvec));
1533   a->vecinuse = col + 1;
1534   PetscCall(MatDenseGetLDA(a->A, &lda));
1535   PetscCall(MatDenseGetArray(a->A, (PetscScalar **)&a->ptrinuse));
1536   PetscCall(VecPlaceArray(a->cvec, PetscSafePointerPlusOffset(a->ptrinuse, (size_t)col * (size_t)lda)));
1537   *v = a->cvec;
1538   PetscFunctionReturn(PETSC_SUCCESS);
1539 }
1540 
1541 PetscErrorCode MatDenseRestoreColumnVec_MPIDense(Mat A, PetscInt col, Vec *v)
1542 {
1543   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1544 
1545   PetscFunctionBegin;
1546   PetscCheck(a->vecinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseGetColumnVec() first");
1547   PetscCheck(a->cvec, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing internal column vector");
1548   VecCheckAssembled(a->cvec);
1549   a->vecinuse = 0;
1550   PetscCall(MatDenseRestoreArray(a->A, (PetscScalar **)&a->ptrinuse));
1551   PetscCall(VecResetArray(a->cvec));
1552   if (v) *v = NULL;
1553   PetscFunctionReturn(PETSC_SUCCESS);
1554 }
1555 
1556 PetscErrorCode MatDenseGetColumnVecRead_MPIDense(Mat A, PetscInt col, Vec *v)
1557 {
1558   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1559   PetscInt      lda;
1560 
1561   PetscFunctionBegin;
1562   PetscCheck(!a->vecinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first");
1563   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
1564   if (!a->cvec) PetscCall(MatDenseCreateColumnVec_Private(A, &a->cvec));
1565   a->vecinuse = col + 1;
1566   PetscCall(MatDenseGetLDA(a->A, &lda));
1567   PetscCall(MatDenseGetArrayRead(a->A, &a->ptrinuse));
1568   PetscCall(VecPlaceArray(a->cvec, PetscSafePointerPlusOffset(a->ptrinuse, (size_t)col * (size_t)lda)));
1569   PetscCall(VecLockReadPush(a->cvec));
1570   *v = a->cvec;
1571   PetscFunctionReturn(PETSC_SUCCESS);
1572 }
1573 
1574 PetscErrorCode MatDenseRestoreColumnVecRead_MPIDense(Mat A, PetscInt col, Vec *v)
1575 {
1576   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1577 
1578   PetscFunctionBegin;
1579   PetscCheck(a->vecinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseGetColumnVec() first");
1580   PetscCheck(a->cvec, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing internal column vector");
1581   VecCheckAssembled(a->cvec);
1582   a->vecinuse = 0;
1583   PetscCall(MatDenseRestoreArrayRead(a->A, &a->ptrinuse));
1584   PetscCall(VecLockReadPop(a->cvec));
1585   PetscCall(VecResetArray(a->cvec));
1586   if (v) *v = NULL;
1587   PetscFunctionReturn(PETSC_SUCCESS);
1588 }
1589 
1590 PetscErrorCode MatDenseGetColumnVecWrite_MPIDense(Mat A, PetscInt col, Vec *v)
1591 {
1592   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1593   PetscInt      lda;
1594 
1595   PetscFunctionBegin;
1596   PetscCheck(!a->vecinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first");
1597   PetscCheck(!a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
1598   if (!a->cvec) PetscCall(MatDenseCreateColumnVec_Private(A, &a->cvec));
1599   a->vecinuse = col + 1;
1600   PetscCall(MatDenseGetLDA(a->A, &lda));
1601   PetscCall(MatDenseGetArrayWrite(a->A, (PetscScalar **)&a->ptrinuse));
1602   PetscCall(VecPlaceArray(a->cvec, PetscSafePointerPlusOffset(a->ptrinuse, (size_t)col * (size_t)lda)));
1603   *v = a->cvec;
1604   PetscFunctionReturn(PETSC_SUCCESS);
1605 }
1606 
1607 PetscErrorCode MatDenseRestoreColumnVecWrite_MPIDense(Mat A, PetscInt col, Vec *v)
1608 {
1609   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1610 
1611   PetscFunctionBegin;
1612   PetscCheck(a->vecinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseGetColumnVec() first");
1613   PetscCheck(a->cvec, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing internal column vector");
1614   VecCheckAssembled(a->cvec);
1615   a->vecinuse = 0;
1616   PetscCall(MatDenseRestoreArrayWrite(a->A, (PetscScalar **)&a->ptrinuse));
1617   PetscCall(VecResetArray(a->cvec));
1618   if (v) *v = NULL;
1619   PetscFunctionReturn(PETSC_SUCCESS);
1620 }
1621 
1622 static PetscErrorCode MatDenseGetSubMatrix_MPIDense(Mat A, PetscInt rbegin, PetscInt rend, PetscInt cbegin, PetscInt cend, Mat *v)
1623 {
1624   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1625   Mat_MPIDense *c;
1626   MPI_Comm      comm;
1627   PetscInt      pbegin, pend;
1628 
1629   PetscFunctionBegin;
1630   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
1631   PetscCheck(!a->vecinuse, comm, PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first");
1632   PetscCheck(!a->matinuse, comm, PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
1633   pbegin = PetscMax(0, PetscMin(A->rmap->rend, rbegin) - A->rmap->rstart);
1634   pend   = PetscMin(A->rmap->n, PetscMax(0, rend - A->rmap->rstart));
1635   if (!a->cmat) {
1636     PetscCall(MatCreate(comm, &a->cmat));
1637     PetscCall(MatSetType(a->cmat, ((PetscObject)A)->type_name));
1638     if (rend - rbegin == A->rmap->N) PetscCall(PetscLayoutReference(A->rmap, &a->cmat->rmap));
1639     else {
1640       PetscCall(PetscLayoutSetLocalSize(a->cmat->rmap, pend - pbegin));
1641       PetscCall(PetscLayoutSetSize(a->cmat->rmap, rend - rbegin));
1642       PetscCall(PetscLayoutSetUp(a->cmat->rmap));
1643     }
1644     PetscCall(PetscLayoutSetSize(a->cmat->cmap, cend - cbegin));
1645     PetscCall(PetscLayoutSetUp(a->cmat->cmap));
1646   } else {
1647     PetscBool same = (PetscBool)(rend - rbegin == a->cmat->rmap->N);
1648     if (same && a->cmat->rmap->N != A->rmap->N) {
1649       same = (PetscBool)(pend - pbegin == a->cmat->rmap->n);
1650       PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &same, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
1651     }
1652     if (!same) {
1653       PetscCall(PetscLayoutDestroy(&a->cmat->rmap));
1654       PetscCall(PetscLayoutCreate(comm, &a->cmat->rmap));
1655       PetscCall(PetscLayoutSetLocalSize(a->cmat->rmap, pend - pbegin));
1656       PetscCall(PetscLayoutSetSize(a->cmat->rmap, rend - rbegin));
1657       PetscCall(PetscLayoutSetUp(a->cmat->rmap));
1658     }
1659     if (cend - cbegin != a->cmat->cmap->N) {
1660       PetscCall(PetscLayoutDestroy(&a->cmat->cmap));
1661       PetscCall(PetscLayoutCreate(comm, &a->cmat->cmap));
1662       PetscCall(PetscLayoutSetSize(a->cmat->cmap, cend - cbegin));
1663       PetscCall(PetscLayoutSetUp(a->cmat->cmap));
1664     }
1665   }
1666   c = (Mat_MPIDense *)a->cmat->data;
1667   PetscCheck(!c->A, comm, PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first");
1668   PetscCall(MatDenseGetSubMatrix(a->A, pbegin, pend, cbegin, cend, &c->A));
1669 
1670   a->cmat->preallocated = PETSC_TRUE;
1671   a->cmat->assembled    = PETSC_TRUE;
1672 #if defined(PETSC_HAVE_DEVICE)
1673   a->cmat->offloadmask = c->A->offloadmask;
1674 #endif
1675   a->matinuse = cbegin + 1;
1676   *v          = a->cmat;
1677   PetscFunctionReturn(PETSC_SUCCESS);
1678 }
1679 
1680 static PetscErrorCode MatDenseRestoreSubMatrix_MPIDense(Mat A, Mat *v)
1681 {
1682   Mat_MPIDense *a = (Mat_MPIDense *)A->data;
1683   Mat_MPIDense *c;
1684 
1685   PetscFunctionBegin;
1686   PetscCheck(a->matinuse, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to call MatDenseGetSubMatrix() first");
1687   PetscCheck(a->cmat, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing internal matrix");
1688   PetscCheck(*v == a->cmat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not the matrix obtained from MatDenseGetSubMatrix()");
1689   a->matinuse = 0;
1690   c           = (Mat_MPIDense *)a->cmat->data;
1691   PetscCall(MatDenseRestoreSubMatrix(a->A, &c->A));
1692   if (v) *v = NULL;
1693 #if defined(PETSC_HAVE_DEVICE)
1694   A->offloadmask = a->A->offloadmask;
1695 #endif
1696   PetscFunctionReturn(PETSC_SUCCESS);
1697 }
1698 
1699 /*MC
1700    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
1701 
1702    Options Database Key:
1703 . -mat_type mpidense - sets the matrix type to `MATMPIDENSE` during a call to `MatSetFromOptions()`
1704 
1705   Level: beginner
1706 
1707 .seealso: [](ch_matrices), `Mat`, `MatCreateDense()`, `MATSEQDENSE`, `MATDENSE`
1708 M*/
1709 PetscErrorCode MatCreate_MPIDense(Mat mat)
1710 {
1711   Mat_MPIDense *a;
1712 
1713   PetscFunctionBegin;
1714   PetscCall(PetscNew(&a));
1715   mat->data   = (void *)a;
1716   mat->ops[0] = MatOps_Values;
1717 
1718   mat->insertmode = NOT_SET_VALUES;
1719 
1720   /* build cache for off array entries formed */
1721   a->donotstash = PETSC_FALSE;
1722 
1723   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)mat), 1, &mat->stash));
1724 
1725   /* stuff used for matrix vector multiply */
1726   a->lvec        = NULL;
1727   a->Mvctx       = NULL;
1728   a->roworiented = PETSC_TRUE;
1729 
1730   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetLDA_C", MatDenseGetLDA_MPIDense));
1731   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseSetLDA_C", MatDenseSetLDA_MPIDense));
1732   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetArray_C", MatDenseGetArray_MPIDense));
1733   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreArray_C", MatDenseRestoreArray_MPIDense));
1734   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetArrayRead_C", MatDenseGetArrayRead_MPIDense));
1735   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreArrayRead_C", MatDenseRestoreArrayRead_MPIDense));
1736   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetArrayWrite_C", MatDenseGetArrayWrite_MPIDense));
1737   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreArrayWrite_C", MatDenseRestoreArrayWrite_MPIDense));
1738   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDensePlaceArray_C", MatDensePlaceArray_MPIDense));
1739   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseResetArray_C", MatDenseResetArray_MPIDense));
1740   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseReplaceArray_C", MatDenseReplaceArray_MPIDense));
1741   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetColumnVec_C", MatDenseGetColumnVec_MPIDense));
1742   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreColumnVec_C", MatDenseRestoreColumnVec_MPIDense));
1743   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetColumnVecRead_C", MatDenseGetColumnVecRead_MPIDense));
1744   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreColumnVecRead_C", MatDenseRestoreColumnVecRead_MPIDense));
1745   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetColumnVecWrite_C", MatDenseGetColumnVecWrite_MPIDense));
1746   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreColumnVecWrite_C", MatDenseRestoreColumnVecWrite_MPIDense));
1747   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetSubMatrix_C", MatDenseGetSubMatrix_MPIDense));
1748   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreSubMatrix_C", MatDenseRestoreSubMatrix_MPIDense));
1749   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpiaij_mpidense_C", MatConvert_MPIAIJ_MPIDense));
1750   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_mpiaij_C", MatConvert_MPIDense_MPIAIJ));
1751 #if defined(PETSC_HAVE_ELEMENTAL)
1752   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_elemental_C", MatConvert_MPIDense_Elemental));
1753 #endif
1754 #if defined(PETSC_HAVE_SCALAPACK)
1755   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_scalapack_C", MatConvert_Dense_ScaLAPACK));
1756 #endif
1757 #if defined(PETSC_HAVE_CUDA)
1758   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_mpidensecuda_C", MatConvert_MPIDense_MPIDenseCUDA));
1759 #endif
1760   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIDenseSetPreallocation_C", MatMPIDenseSetPreallocation_MPIDense));
1761   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaij_mpidense_C", MatProductSetFromOptions_MPIAIJ_MPIDense));
1762   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidense_mpiaij_C", MatProductSetFromOptions_MPIDense_MPIAIJ));
1763 #if defined(PETSC_HAVE_CUDA)
1764   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaijcusparse_mpidense_C", MatProductSetFromOptions_MPIAIJ_MPIDense));
1765   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidense_mpiaijcusparse_C", MatProductSetFromOptions_MPIDense_MPIAIJ));
1766 #endif
1767 #if defined(PETSC_HAVE_HIP)
1768   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpidense_mpidensehip_C", MatConvert_MPIDense_MPIDenseHIP));
1769   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpiaijhipsparse_mpidense_C", MatProductSetFromOptions_MPIAIJ_MPIDense));
1770   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_mpidense_mpiaijhipsparse_C", MatProductSetFromOptions_MPIDense_MPIAIJ));
1771 #endif
1772   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseGetColumn_C", MatDenseGetColumn_MPIDense));
1773   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDenseRestoreColumn_C", MatDenseRestoreColumn_MPIDense));
1774   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMultAddColumnRange_C", MatMultAddColumnRange_MPIDense));
1775   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMultHermitianTransposeColumnRange_C", MatMultHermitianTransposeColumnRange_MPIDense));
1776   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMultHermitianTransposeAddColumnRange_C", MatMultHermitianTransposeAddColumnRange_MPIDense));
1777   PetscCall(PetscObjectChangeTypeName((PetscObject)mat, MATMPIDENSE));
1778   PetscFunctionReturn(PETSC_SUCCESS);
1779 }
1780 
1781 /*MC
1782    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1783 
1784    This matrix type is identical to `MATSEQDENSE` when constructed with a single process communicator,
1785    and `MATMPIDENSE` otherwise.
1786 
1787    Options Database Key:
1788 . -mat_type dense - sets the matrix type to `MATDENSE` during a call to `MatSetFromOptions()`
1789 
1790   Level: beginner
1791 
1792 .seealso: [](ch_matrices), `Mat`, `MATSEQDENSE`, `MATMPIDENSE`, `MATDENSECUDA`, `MATDENSEHIP`
1793 M*/
1794 
1795 /*@
1796   MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1797 
1798   Collective
1799 
1800   Input Parameters:
1801 + B    - the matrix
1802 - data - optional location of matrix data.  Set to `NULL` for PETSc
1803          to control all matrix memory allocation.
1804 
1805   Level: intermediate
1806 
1807   Notes:
1808   The dense format is fully compatible with standard Fortran
1809   storage by columns.
1810 
1811   The data input variable is intended primarily for Fortran programmers
1812   who wish to allocate their own matrix memory space.  Most users should
1813   set `data` to `NULL`.
1814 
1815 .seealso: [](ch_matrices), `Mat`, `MATMPIDENSE`, `MatCreate()`, `MatCreateSeqDense()`, `MatSetValues()`
1816 @*/
1817 PetscErrorCode MatMPIDenseSetPreallocation(Mat B, PetscScalar *data)
1818 {
1819   PetscFunctionBegin;
1820   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
1821   PetscTryMethod(B, "MatMPIDenseSetPreallocation_C", (Mat, PetscScalar *), (B, data));
1822   PetscFunctionReturn(PETSC_SUCCESS);
1823 }
1824 
1825 /*@
1826   MatDensePlaceArray - Allows one to replace the array in a `MATDENSE` matrix with an
1827   array provided by the user. This is useful to avoid copying an array
1828   into a matrix
1829 
1830   Not Collective
1831 
1832   Input Parameters:
1833 + mat   - the matrix
1834 - array - the array in column major order
1835 
1836   Level: developer
1837 
1838   Note:
1839   You can return to the original array with a call to `MatDenseResetArray()`. The user is responsible for freeing this array; it will not be
1840   freed when the matrix is destroyed.
1841 
1842 .seealso: [](ch_matrices), `Mat`, `MATDENSE`, `MatDenseGetArray()`, `MatDenseResetArray()`, `VecPlaceArray()`, `VecGetArray()`, `VecRestoreArray()`, `VecReplaceArray()`, `VecResetArray()`,
1843           `MatDenseReplaceArray()`
1844 @*/
1845 PetscErrorCode MatDensePlaceArray(Mat mat, const PetscScalar *array)
1846 {
1847   PetscFunctionBegin;
1848   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
1849   PetscUseMethod(mat, "MatDensePlaceArray_C", (Mat, const PetscScalar *), (mat, array));
1850   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
1851 #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
1852   mat->offloadmask = PETSC_OFFLOAD_CPU;
1853 #endif
1854   PetscFunctionReturn(PETSC_SUCCESS);
1855 }
1856 
1857 /*@
1858   MatDenseResetArray - Resets the matrix array to that it previously had before the call to `MatDensePlaceArray()`
1859 
1860   Not Collective
1861 
1862   Input Parameter:
1863 . mat - the matrix
1864 
1865   Level: developer
1866 
1867   Note:
1868   You can only call this after a call to `MatDensePlaceArray()`
1869 
1870 .seealso: [](ch_matrices), `Mat`, `MATDENSE`, `MatDenseGetArray()`, `MatDensePlaceArray()`, `VecPlaceArray()`, `VecGetArray()`, `VecRestoreArray()`, `VecReplaceArray()`, `VecResetArray()`
1871 @*/
1872 PetscErrorCode MatDenseResetArray(Mat mat)
1873 {
1874   PetscFunctionBegin;
1875   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
1876   PetscUseMethod(mat, "MatDenseResetArray_C", (Mat), (mat));
1877   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
1878   PetscFunctionReturn(PETSC_SUCCESS);
1879 }
1880 
1881 /*@
1882   MatDenseReplaceArray - Allows one to replace the array in a dense matrix with an
1883   array provided by the user. This is useful to avoid copying an array
1884   into a matrix
1885 
1886   Not Collective
1887 
1888   Input Parameters:
1889 + mat   - the matrix
1890 - array - the array in column major order
1891 
1892   Level: developer
1893 
1894   Note:
1895   The memory passed in MUST be obtained with `PetscMalloc()` and CANNOT be
1896   freed by the user. It will be freed when the matrix is destroyed.
1897 
1898 .seealso: [](ch_matrices), `Mat`, `MatDensePlaceArray()`, `MatDenseGetArray()`, `VecReplaceArray()`
1899 @*/
1900 PetscErrorCode MatDenseReplaceArray(Mat mat, const PetscScalar *array)
1901 {
1902   PetscFunctionBegin;
1903   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
1904   PetscUseMethod(mat, "MatDenseReplaceArray_C", (Mat, const PetscScalar *), (mat, array));
1905   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
1906 #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
1907   mat->offloadmask = PETSC_OFFLOAD_CPU;
1908 #endif
1909   PetscFunctionReturn(PETSC_SUCCESS);
1910 }
1911 
1912 /*@
1913   MatCreateDense - Creates a matrix in `MATDENSE` format.
1914 
1915   Collective
1916 
1917   Input Parameters:
1918 + comm - MPI communicator
1919 . m    - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
1920 . n    - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
1921 . M    - number of global rows (or `PETSC_DECIDE` to have calculated if `m` is given)
1922 . N    - number of global columns (or `PETSC_DECIDE` to have calculated if `n` is given)
1923 - data - optional location of matrix data.  Set data to `NULL` (`PETSC_NULL_SCALAR` for Fortran users) for PETSc
1924    to control all matrix memory allocation.
1925 
1926   Output Parameter:
1927 . A - the matrix
1928 
1929   Level: intermediate
1930 
1931   Notes:
1932   The dense format is fully compatible with standard Fortran
1933   storage by columns.
1934 
1935   Although local portions of the matrix are stored in column-major
1936   order, the matrix is partitioned across MPI ranks by row.
1937 
1938   The data input variable is intended primarily for Fortran programmers
1939   who wish to allocate their own matrix memory space.  Most users should
1940   set `data` to `NULL` (`PETSC_NULL_SCALAR` for Fortran users).
1941 
1942   The user MUST specify either the local or global matrix dimensions
1943   (possibly both).
1944 
1945 .seealso: [](ch_matrices), `Mat`, `MATDENSE`, `MatCreate()`, `MatCreateSeqDense()`, `MatSetValues()`
1946 @*/
1947 PetscErrorCode MatCreateDense(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar *data, Mat *A)
1948 {
1949   PetscFunctionBegin;
1950   PetscCall(MatCreate(comm, A));
1951   PetscCall(MatSetSizes(*A, m, n, M, N));
1952   PetscCall(MatSetType(*A, MATDENSE));
1953   PetscCall(MatSeqDenseSetPreallocation(*A, data));
1954   PetscCall(MatMPIDenseSetPreallocation(*A, data));
1955   PetscFunctionReturn(PETSC_SUCCESS);
1956 }
1957 
1958 static PetscErrorCode MatDuplicate_MPIDense(Mat A, MatDuplicateOption cpvalues, Mat *newmat)
1959 {
1960   Mat           mat;
1961   Mat_MPIDense *a, *oldmat = (Mat_MPIDense *)A->data;
1962 
1963   PetscFunctionBegin;
1964   *newmat = NULL;
1965   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat));
1966   PetscCall(MatSetSizes(mat, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
1967   PetscCall(MatSetType(mat, ((PetscObject)A)->type_name));
1968   a = (Mat_MPIDense *)mat->data;
1969 
1970   mat->factortype   = A->factortype;
1971   mat->assembled    = PETSC_TRUE;
1972   mat->preallocated = PETSC_TRUE;
1973 
1974   mat->insertmode = NOT_SET_VALUES;
1975   a->donotstash   = oldmat->donotstash;
1976 
1977   PetscCall(PetscLayoutReference(A->rmap, &mat->rmap));
1978   PetscCall(PetscLayoutReference(A->cmap, &mat->cmap));
1979 
1980   PetscCall(MatDuplicate(oldmat->A, cpvalues, &a->A));
1981 
1982   *newmat = mat;
1983   PetscFunctionReturn(PETSC_SUCCESS);
1984 }
1985 
1986 static PetscErrorCode MatLoad_MPIDense(Mat newMat, PetscViewer viewer)
1987 {
1988   PetscBool isbinary;
1989 #if defined(PETSC_HAVE_HDF5)
1990   PetscBool ishdf5;
1991 #endif
1992 
1993   PetscFunctionBegin;
1994   PetscValidHeaderSpecific(newMat, MAT_CLASSID, 1);
1995   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1996   /* force binary viewer to load .info file if it has not yet done so */
1997   PetscCall(PetscViewerSetUp(viewer));
1998   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1999 #if defined(PETSC_HAVE_HDF5)
2000   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2001 #endif
2002   if (isbinary) {
2003     PetscCall(MatLoad_Dense_Binary(newMat, viewer));
2004 #if defined(PETSC_HAVE_HDF5)
2005   } else if (ishdf5) {
2006     PetscCall(MatLoad_Dense_HDF5(newMat, viewer));
2007 #endif
2008   } else SETERRQ(PetscObjectComm((PetscObject)newMat), PETSC_ERR_SUP, "Viewer type %s not yet supported for reading %s matrices", ((PetscObject)viewer)->type_name, ((PetscObject)newMat)->type_name);
2009   PetscFunctionReturn(PETSC_SUCCESS);
2010 }
2011 
2012 static PetscErrorCode MatEqual_MPIDense(Mat A, Mat B, PetscBool *flag)
2013 {
2014   Mat_MPIDense *matB = (Mat_MPIDense *)B->data, *matA = (Mat_MPIDense *)A->data;
2015   Mat           a, b;
2016 
2017   PetscFunctionBegin;
2018   a = matA->A;
2019   b = matB->A;
2020   PetscCall(MatEqual(a, b, flag));
2021   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, flag, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
2022   PetscFunctionReturn(PETSC_SUCCESS);
2023 }
2024 
2025 static PetscErrorCode MatDestroy_MatTransMatMult_MPIDense_MPIDense(void *data)
2026 {
2027   Mat_TransMatMultDense *atb = (Mat_TransMatMultDense *)data;
2028 
2029   PetscFunctionBegin;
2030   PetscCall(PetscFree2(atb->sendbuf, atb->recvcounts));
2031   PetscCall(MatDestroy(&atb->atb));
2032   PetscCall(PetscFree(atb));
2033   PetscFunctionReturn(PETSC_SUCCESS);
2034 }
2035 
2036 static PetscErrorCode MatDestroy_MatMatTransMult_MPIDense_MPIDense(void *data)
2037 {
2038   Mat_MatTransMultDense *abt = (Mat_MatTransMultDense *)data;
2039 
2040   PetscFunctionBegin;
2041   PetscCall(PetscFree2(abt->buf[0], abt->buf[1]));
2042   PetscCall(PetscFree2(abt->recvcounts, abt->recvdispls));
2043   PetscCall(PetscFree(abt));
2044   PetscFunctionReturn(PETSC_SUCCESS);
2045 }
2046 
2047 static PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat A, Mat B, Mat C)
2048 {
2049   Mat_MPIDense          *a = (Mat_MPIDense *)A->data, *b = (Mat_MPIDense *)B->data, *c = (Mat_MPIDense *)C->data;
2050   Mat_TransMatMultDense *atb;
2051   MPI_Comm               comm;
2052   PetscMPIInt            size, *recvcounts;
2053   PetscScalar           *carray, *sendbuf;
2054   const PetscScalar     *atbarray;
2055   PetscInt               i, cN = C->cmap->N, proc, k, j, lda;
2056   const PetscInt        *ranges;
2057 
2058   PetscFunctionBegin;
2059   MatCheckProduct(C, 3);
2060   PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
2061   atb        = (Mat_TransMatMultDense *)C->product->data;
2062   recvcounts = atb->recvcounts;
2063   sendbuf    = atb->sendbuf;
2064 
2065   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
2066   PetscCallMPI(MPI_Comm_size(comm, &size));
2067 
2068   /* compute atbarray = aseq^T * bseq */
2069   PetscCall(MatTransposeMatMult(a->A, b->A, atb->atb ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX, PETSC_DETERMINE, &atb->atb));
2070 
2071   PetscCall(MatGetOwnershipRanges(C, &ranges));
2072 
2073   /* arrange atbarray into sendbuf */
2074   PetscCall(MatDenseGetArrayRead(atb->atb, &atbarray));
2075   PetscCall(MatDenseGetLDA(atb->atb, &lda));
2076   for (proc = 0, k = 0; proc < size; proc++) {
2077     for (j = 0; j < cN; j++) {
2078       for (i = ranges[proc]; i < ranges[proc + 1]; i++) sendbuf[k++] = atbarray[i + j * lda];
2079     }
2080   }
2081   PetscCall(MatDenseRestoreArrayRead(atb->atb, &atbarray));
2082 
2083   /* sum all atbarray to local values of C */
2084   PetscCall(MatDenseGetArrayWrite(c->A, &carray));
2085   PetscCallMPI(MPI_Reduce_scatter(sendbuf, carray, recvcounts, MPIU_SCALAR, MPIU_SUM, comm));
2086   PetscCall(MatDenseRestoreArrayWrite(c->A, &carray));
2087   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
2088   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2089   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
2090   PetscFunctionReturn(PETSC_SUCCESS);
2091 }
2092 
2093 static PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat A, Mat B, PetscReal fill, Mat C)
2094 {
2095   MPI_Comm               comm;
2096   PetscMPIInt            size;
2097   PetscInt               cm = A->cmap->n, cM, cN = B->cmap->N;
2098   Mat_TransMatMultDense *atb;
2099   PetscBool              cisdense = PETSC_FALSE;
2100   PetscInt               i;
2101   const PetscInt        *ranges;
2102 
2103   PetscFunctionBegin;
2104   MatCheckProduct(C, 4);
2105   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
2106   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
2107   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) {
2108     SETERRQ(comm, PETSC_ERR_ARG_SIZ, "Matrix local dimensions are incompatible, A (%" PetscInt_FMT ", %" PetscInt_FMT ") != B (%" PetscInt_FMT ",%" PetscInt_FMT ")", A->rmap->rstart, A->rmap->rend, B->rmap->rstart, B->rmap->rend);
2109   }
2110 
2111   /* create matrix product C */
2112   PetscCall(MatSetSizes(C, cm, B->cmap->n, A->cmap->N, B->cmap->N));
2113 #if defined(PETSC_HAVE_CUDA)
2114   PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATMPIDENSE, MATMPIDENSECUDA, ""));
2115 #endif
2116 #if defined(PETSC_HAVE_HIP)
2117   PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATMPIDENSE, MATMPIDENSEHIP, ""));
2118 #endif
2119   if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
2120   PetscCall(MatSetUp(C));
2121 
2122   /* create data structure for reuse C */
2123   PetscCallMPI(MPI_Comm_size(comm, &size));
2124   PetscCall(PetscNew(&atb));
2125   cM = C->rmap->N;
2126   PetscCall(PetscMalloc2(cM * cN, &atb->sendbuf, size, &atb->recvcounts));
2127   PetscCall(MatGetOwnershipRanges(C, &ranges));
2128   for (i = 0; i < size; i++) atb->recvcounts[i] = (ranges[i + 1] - ranges[i]) * cN;
2129 
2130   C->product->data    = atb;
2131   C->product->destroy = MatDestroy_MatTransMatMult_MPIDense_MPIDense;
2132   PetscFunctionReturn(PETSC_SUCCESS);
2133 }
2134 
2135 static PetscErrorCode MatMatTransposeMultSymbolic_MPIDense_MPIDense(Mat A, Mat B, PetscReal fill, Mat C)
2136 {
2137   MPI_Comm               comm;
2138   PetscMPIInt            i, size;
2139   PetscInt               maxRows, bufsiz;
2140   PetscMPIInt            tag;
2141   PetscInt               alg;
2142   Mat_MatTransMultDense *abt;
2143   Mat_Product           *product = C->product;
2144   PetscBool              flg;
2145 
2146   PetscFunctionBegin;
2147   MatCheckProduct(C, 4);
2148   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
2149   /* check local size of A and B */
2150   PetscCheck(A->cmap->n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrix local column dimensions are incompatible, A (%" PetscInt_FMT ") != B (%" PetscInt_FMT ")", A->cmap->n, B->cmap->n);
2151 
2152   PetscCall(PetscStrcmp(product->alg, "allgatherv", &flg));
2153   alg = flg ? 0 : 1;
2154 
2155   /* setup matrix product C */
2156   PetscCall(MatSetSizes(C, A->rmap->n, B->rmap->n, A->rmap->N, B->rmap->N));
2157   PetscCall(MatSetType(C, MATMPIDENSE));
2158   PetscCall(MatSetUp(C));
2159   PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag));
2160 
2161   /* create data structure for reuse C */
2162   PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
2163   PetscCallMPI(MPI_Comm_size(comm, &size));
2164   PetscCall(PetscNew(&abt));
2165   abt->tag = tag;
2166   abt->alg = alg;
2167   switch (alg) {
2168   case 1: /* alg: "cyclic" */
2169     for (maxRows = 0, i = 0; i < size; i++) maxRows = PetscMax(maxRows, (B->rmap->range[i + 1] - B->rmap->range[i]));
2170     bufsiz = A->cmap->N * maxRows;
2171     PetscCall(PetscMalloc2(bufsiz, &abt->buf[0], bufsiz, &abt->buf[1]));
2172     break;
2173   default: /* alg: "allgatherv" */
2174     PetscCall(PetscMalloc2(B->rmap->n * B->cmap->N, &abt->buf[0], B->rmap->N * B->cmap->N, &abt->buf[1]));
2175     PetscCall(PetscMalloc2(size, &abt->recvcounts, size + 1, &abt->recvdispls));
2176     for (i = 0; i <= size; i++) abt->recvdispls[i] = B->rmap->range[i] * A->cmap->N;
2177     for (i = 0; i < size; i++) abt->recvcounts[i] = abt->recvdispls[i + 1] - abt->recvdispls[i];
2178     break;
2179   }
2180 
2181   C->product->data                = abt;
2182   C->product->destroy             = MatDestroy_MatMatTransMult_MPIDense_MPIDense;
2183   C->ops->mattransposemultnumeric = MatMatTransposeMultNumeric_MPIDense_MPIDense;
2184   PetscFunctionReturn(PETSC_SUCCESS);
2185 }
2186 
2187 static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense_Cyclic(Mat A, Mat B, Mat C)
2188 {
2189   Mat_MPIDense          *a = (Mat_MPIDense *)A->data, *b = (Mat_MPIDense *)B->data, *c = (Mat_MPIDense *)C->data;
2190   Mat_MatTransMultDense *abt;
2191   MPI_Comm               comm;
2192   PetscMPIInt            rank, size, sendsiz, recvsiz, sendto, recvfrom, recvisfrom;
2193   PetscScalar           *sendbuf, *recvbuf = NULL, *cv;
2194   PetscInt               i, cK             = A->cmap->N, k, j, bn;
2195   PetscScalar            _DOne = 1.0, _DZero = 0.0;
2196   const PetscScalar     *av, *bv;
2197   PetscBLASInt           cm, cn, ck, alda, blda = 0, clda;
2198   MPI_Request            reqs[2];
2199   const PetscInt        *ranges;
2200 
2201   PetscFunctionBegin;
2202   MatCheckProduct(C, 3);
2203   PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
2204   abt = (Mat_MatTransMultDense *)C->product->data;
2205   PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
2206   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2207   PetscCallMPI(MPI_Comm_size(comm, &size));
2208   PetscCall(MatDenseGetArrayRead(a->A, &av));
2209   PetscCall(MatDenseGetArrayRead(b->A, &bv));
2210   PetscCall(MatDenseGetArrayWrite(c->A, &cv));
2211   PetscCall(MatDenseGetLDA(a->A, &i));
2212   PetscCall(PetscBLASIntCast(i, &alda));
2213   PetscCall(MatDenseGetLDA(b->A, &i));
2214   PetscCall(PetscBLASIntCast(i, &blda));
2215   PetscCall(MatDenseGetLDA(c->A, &i));
2216   PetscCall(PetscBLASIntCast(i, &clda));
2217   PetscCall(MatGetOwnershipRanges(B, &ranges));
2218   bn = B->rmap->n;
2219   if (blda == bn) {
2220     sendbuf = (PetscScalar *)bv;
2221   } else {
2222     sendbuf = abt->buf[0];
2223     for (k = 0, i = 0; i < cK; i++) {
2224       for (j = 0; j < bn; j++, k++) sendbuf[k] = bv[i * blda + j];
2225     }
2226   }
2227   if (size > 1) {
2228     sendto   = (rank + size - 1) % size;
2229     recvfrom = (rank + size + 1) % size;
2230   } else {
2231     sendto = recvfrom = 0;
2232   }
2233   PetscCall(PetscBLASIntCast(cK, &ck));
2234   PetscCall(PetscBLASIntCast(c->A->rmap->n, &cm));
2235   recvisfrom = rank;
2236   for (i = 0; i < size; i++) {
2237     /* we have finished receiving in sending, bufs can be read/modified */
2238     PetscInt nextrecvisfrom = (recvisfrom + 1) % size; /* which process the next recvbuf will originate on */
2239     PetscInt nextbn         = ranges[nextrecvisfrom + 1] - ranges[nextrecvisfrom];
2240 
2241     if (nextrecvisfrom != rank) {
2242       /* start the cyclic sends from sendbuf, to recvbuf (which will switch to sendbuf) */
2243       sendsiz = cK * bn;
2244       recvsiz = cK * nextbn;
2245       recvbuf = (i & 1) ? abt->buf[0] : abt->buf[1];
2246       PetscCallMPI(MPI_Isend(sendbuf, sendsiz, MPIU_SCALAR, sendto, abt->tag, comm, &reqs[0]));
2247       PetscCallMPI(MPI_Irecv(recvbuf, recvsiz, MPIU_SCALAR, recvfrom, abt->tag, comm, &reqs[1]));
2248     }
2249 
2250     /* local aseq * sendbuf^T */
2251     PetscCall(PetscBLASIntCast(ranges[recvisfrom + 1] - ranges[recvisfrom], &cn));
2252     if (cm && cn && ck) PetscCallBLAS("BLASgemm", BLASgemm_("N", "T", &cm, &cn, &ck, &_DOne, av, &alda, sendbuf, &cn, &_DZero, cv + clda * ranges[recvisfrom], &clda));
2253 
2254     if (nextrecvisfrom != rank) {
2255       /* wait for the sends and receives to complete, swap sendbuf and recvbuf */
2256       PetscCallMPI(MPI_Waitall(2, reqs, MPI_STATUSES_IGNORE));
2257     }
2258     bn         = nextbn;
2259     recvisfrom = nextrecvisfrom;
2260     sendbuf    = recvbuf;
2261   }
2262   PetscCall(MatDenseRestoreArrayRead(a->A, &av));
2263   PetscCall(MatDenseRestoreArrayRead(b->A, &bv));
2264   PetscCall(MatDenseRestoreArrayWrite(c->A, &cv));
2265   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
2266   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2267   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
2268   PetscFunctionReturn(PETSC_SUCCESS);
2269 }
2270 
2271 static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense_Allgatherv(Mat A, Mat B, Mat C)
2272 {
2273   Mat_MPIDense          *a = (Mat_MPIDense *)A->data, *b = (Mat_MPIDense *)B->data, *c = (Mat_MPIDense *)C->data;
2274   Mat_MatTransMultDense *abt;
2275   MPI_Comm               comm;
2276   PetscMPIInt            size;
2277   PetscScalar           *cv, *sendbuf, *recvbuf;
2278   const PetscScalar     *av, *bv;
2279   PetscInt               blda, i, cK = A->cmap->N, k, j, bn;
2280   PetscScalar            _DOne = 1.0, _DZero = 0.0;
2281   PetscBLASInt           cm, cn, ck, alda, clda;
2282 
2283   PetscFunctionBegin;
2284   MatCheckProduct(C, 3);
2285   PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
2286   abt = (Mat_MatTransMultDense *)C->product->data;
2287   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
2288   PetscCallMPI(MPI_Comm_size(comm, &size));
2289   PetscCall(MatDenseGetArrayRead(a->A, &av));
2290   PetscCall(MatDenseGetArrayRead(b->A, &bv));
2291   PetscCall(MatDenseGetArrayWrite(c->A, &cv));
2292   PetscCall(MatDenseGetLDA(a->A, &i));
2293   PetscCall(PetscBLASIntCast(i, &alda));
2294   PetscCall(MatDenseGetLDA(b->A, &blda));
2295   PetscCall(MatDenseGetLDA(c->A, &i));
2296   PetscCall(PetscBLASIntCast(i, &clda));
2297   /* copy transpose of B into buf[0] */
2298   bn      = B->rmap->n;
2299   sendbuf = abt->buf[0];
2300   recvbuf = abt->buf[1];
2301   for (k = 0, j = 0; j < bn; j++) {
2302     for (i = 0; i < cK; i++, k++) sendbuf[k] = bv[i * blda + j];
2303   }
2304   PetscCall(MatDenseRestoreArrayRead(b->A, &bv));
2305   PetscCallMPI(MPI_Allgatherv(sendbuf, bn * cK, MPIU_SCALAR, recvbuf, abt->recvcounts, abt->recvdispls, MPIU_SCALAR, comm));
2306   PetscCall(PetscBLASIntCast(cK, &ck));
2307   PetscCall(PetscBLASIntCast(c->A->rmap->n, &cm));
2308   PetscCall(PetscBLASIntCast(c->A->cmap->n, &cn));
2309   if (cm && cn && ck) PetscCallBLAS("BLASgemm", BLASgemm_("N", "N", &cm, &cn, &ck, &_DOne, av, &alda, recvbuf, &ck, &_DZero, cv, &clda));
2310   PetscCall(MatDenseRestoreArrayRead(a->A, &av));
2311   PetscCall(MatDenseRestoreArrayRead(b->A, &bv));
2312   PetscCall(MatDenseRestoreArrayWrite(c->A, &cv));
2313   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
2314   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2315   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
2316   PetscFunctionReturn(PETSC_SUCCESS);
2317 }
2318 
2319 static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense(Mat A, Mat B, Mat C)
2320 {
2321   Mat_MatTransMultDense *abt;
2322 
2323   PetscFunctionBegin;
2324   MatCheckProduct(C, 3);
2325   PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
2326   abt = (Mat_MatTransMultDense *)C->product->data;
2327   switch (abt->alg) {
2328   case 1:
2329     PetscCall(MatMatTransposeMultNumeric_MPIDense_MPIDense_Cyclic(A, B, C));
2330     break;
2331   default:
2332     PetscCall(MatMatTransposeMultNumeric_MPIDense_MPIDense_Allgatherv(A, B, C));
2333     break;
2334   }
2335   PetscFunctionReturn(PETSC_SUCCESS);
2336 }
2337 
2338 static PetscErrorCode MatDestroy_MatMatMult_MPIDense_MPIDense(void *data)
2339 {
2340   Mat_MatMultDense *ab = (Mat_MatMultDense *)data;
2341 
2342   PetscFunctionBegin;
2343   PetscCall(MatDestroy(&ab->Ce));
2344   PetscCall(MatDestroy(&ab->Ae));
2345   PetscCall(MatDestroy(&ab->Be));
2346   PetscCall(PetscFree(ab));
2347   PetscFunctionReturn(PETSC_SUCCESS);
2348 }
2349 
2350 static PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat A, Mat B, Mat C)
2351 {
2352   Mat_MatMultDense *ab;
2353   Mat_MPIDense     *mdn = (Mat_MPIDense *)A->data;
2354 
2355   PetscFunctionBegin;
2356   MatCheckProduct(C, 3);
2357   PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing product data");
2358   ab = (Mat_MatMultDense *)C->product->data;
2359   if (ab->Ae && ab->Ce) {
2360 #if PetscDefined(HAVE_ELEMENTAL)
2361     PetscCall(MatConvert_MPIDense_Elemental(A, MATELEMENTAL, MAT_REUSE_MATRIX, &ab->Ae));
2362     PetscCall(MatConvert_MPIDense_Elemental(B, MATELEMENTAL, MAT_REUSE_MATRIX, &ab->Be));
2363     PetscCall(MatMatMultNumeric_Elemental(ab->Ae, ab->Be, ab->Ce));
2364     PetscCall(MatConvert(ab->Ce, MATMPIDENSE, MAT_REUSE_MATRIX, &C));
2365 #else
2366     SETERRQ(PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "PETSC_HAVE_ELEMENTAL not defined");
2367 #endif
2368   } else {
2369     const PetscScalar *read;
2370     PetscScalar       *write;
2371     PetscInt           lda;
2372 
2373     PetscCall(MatDenseGetLDA(B, &lda));
2374     PetscCall(MatDenseGetArrayRead(B, &read));
2375     PetscCall(MatDenseGetArrayWrite(ab->Be, &write));
2376     if (!mdn->Mvctx) PetscCall(MatSetUpMultiply_MPIDense(A)); /* cannot be done during the symbolic phase because of possible calls to MatProductReplaceMats() */
2377     for (PetscInt i = 0; i < C->cmap->N; ++i) {
2378       PetscCall(PetscSFBcastBegin(mdn->Mvctx, MPIU_SCALAR, read + i * lda, write + i * ab->Be->rmap->n, MPI_REPLACE));
2379       PetscCall(PetscSFBcastEnd(mdn->Mvctx, MPIU_SCALAR, read + i * lda, write + i * ab->Be->rmap->n, MPI_REPLACE));
2380     }
2381     PetscCall(MatDenseRestoreArrayWrite(ab->Be, &write));
2382     PetscCall(MatDenseRestoreArrayRead(B, &read));
2383     PetscCall(MatMatMultNumeric_SeqDense_SeqDense(((Mat_MPIDense *)A->data)->A, ab->Be, ((Mat_MPIDense *)C->data)->A));
2384   }
2385   PetscFunctionReturn(PETSC_SUCCESS);
2386 }
2387 
2388 static PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A, Mat B, PetscReal fill, Mat C)
2389 {
2390   Mat_Product      *product = C->product;
2391   PetscInt          alg;
2392   Mat_MatMultDense *ab;
2393   PetscBool         flg;
2394 
2395   PetscFunctionBegin;
2396   MatCheckProduct(C, 4);
2397   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
2398   /* check local size of A and B */
2399   PetscCheck(A->cmap->rstart == B->rmap->rstart && A->cmap->rend == B->rmap->rend, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Matrix local dimensions are incompatible, A (%" PetscInt_FMT ", %" PetscInt_FMT ") != B (%" PetscInt_FMT ", %" PetscInt_FMT ")",
2400              A->rmap->rstart, A->rmap->rend, B->rmap->rstart, B->rmap->rend);
2401 
2402   PetscCall(PetscStrcmp(product->alg, "petsc", &flg));
2403   alg = flg ? 0 : 1;
2404 
2405   /* setup C */
2406   PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N));
2407   PetscCall(MatSetType(C, MATMPIDENSE));
2408   PetscCall(MatSetUp(C));
2409 
2410   /* create data structure for reuse Cdense */
2411   PetscCall(PetscNew(&ab));
2412 
2413   switch (alg) {
2414   case 1: /* alg: "elemental" */
2415 #if PetscDefined(HAVE_ELEMENTAL)
2416     /* create elemental matrices Ae and Be */
2417     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &ab->Ae));
2418     PetscCall(MatSetSizes(ab->Ae, PETSC_DECIDE, PETSC_DECIDE, A->rmap->N, A->cmap->N));
2419     PetscCall(MatSetType(ab->Ae, MATELEMENTAL));
2420     PetscCall(MatSetUp(ab->Ae));
2421     PetscCall(MatSetOption(ab->Ae, MAT_ROW_ORIENTED, PETSC_FALSE));
2422 
2423     PetscCall(MatCreate(PetscObjectComm((PetscObject)B), &ab->Be));
2424     PetscCall(MatSetSizes(ab->Be, PETSC_DECIDE, PETSC_DECIDE, B->rmap->N, B->cmap->N));
2425     PetscCall(MatSetType(ab->Be, MATELEMENTAL));
2426     PetscCall(MatSetUp(ab->Be));
2427     PetscCall(MatSetOption(ab->Be, MAT_ROW_ORIENTED, PETSC_FALSE));
2428 
2429     /* compute symbolic Ce = Ae*Be */
2430     PetscCall(MatCreate(PetscObjectComm((PetscObject)C), &ab->Ce));
2431     PetscCall(MatMatMultSymbolic_Elemental(ab->Ae, ab->Be, fill, ab->Ce));
2432 #else
2433     SETERRQ(PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "PETSC_HAVE_ELEMENTAL not defined");
2434 #endif
2435     break;
2436   default: /* alg: "petsc" */
2437     ab->Ae = NULL;
2438     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, A->cmap->N, B->cmap->N, NULL, &ab->Be));
2439     ab->Ce = NULL;
2440     break;
2441   }
2442 
2443   C->product->data       = ab;
2444   C->product->destroy    = MatDestroy_MatMatMult_MPIDense_MPIDense;
2445   C->ops->matmultnumeric = MatMatMultNumeric_MPIDense_MPIDense;
2446   PetscFunctionReturn(PETSC_SUCCESS);
2447 }
2448 
2449 static PetscErrorCode MatProductSetFromOptions_MPIDense_AB(Mat C)
2450 {
2451   Mat_Product *product     = C->product;
2452   const char  *algTypes[2] = {"petsc", "elemental"};
2453   PetscInt     alg, nalg = PetscDefined(HAVE_ELEMENTAL) ? 2 : 1;
2454   PetscBool    flg = PETSC_FALSE;
2455 
2456   PetscFunctionBegin;
2457   /* Set default algorithm */
2458   alg = 0; /* default is petsc */
2459   PetscCall(PetscStrcmp(product->alg, "default", &flg));
2460   if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2461 
2462   /* Get runtime option */
2463   PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatProduct_AB", "Mat");
2464   PetscCall(PetscOptionsEList("-mat_product_algorithm", "Algorithmic approach", "MatProduct_AB", algTypes, nalg, algTypes[alg], &alg, &flg));
2465   PetscOptionsEnd();
2466   if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2467 
2468   C->ops->matmultsymbolic = MatMatMultSymbolic_MPIDense_MPIDense;
2469   C->ops->productsymbolic = MatProductSymbolic_AB;
2470   PetscFunctionReturn(PETSC_SUCCESS);
2471 }
2472 
2473 static PetscErrorCode MatProductSetFromOptions_MPIDense_AtB(Mat C)
2474 {
2475   Mat_Product *product = C->product;
2476   Mat          A = product->A, B = product->B;
2477 
2478   PetscFunctionBegin;
2479   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, (%" PetscInt_FMT ", %" PetscInt_FMT ") != (%" PetscInt_FMT ",%" PetscInt_FMT ")",
2480              A->rmap->rstart, A->rmap->rend, B->rmap->rstart, B->rmap->rend);
2481   C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_MPIDense_MPIDense;
2482   C->ops->productsymbolic          = MatProductSymbolic_AtB;
2483   PetscFunctionReturn(PETSC_SUCCESS);
2484 }
2485 
2486 static PetscErrorCode MatProductSetFromOptions_MPIDense_ABt(Mat C)
2487 {
2488   Mat_Product *product     = C->product;
2489   const char  *algTypes[2] = {"allgatherv", "cyclic"};
2490   PetscInt     alg, nalg = 2;
2491   PetscBool    flg = PETSC_FALSE;
2492 
2493   PetscFunctionBegin;
2494   /* Set default algorithm */
2495   alg = 0; /* default is allgatherv */
2496   PetscCall(PetscStrcmp(product->alg, "default", &flg));
2497   if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2498 
2499   /* Get runtime option */
2500   if (product->api_user) {
2501     PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatMatTransposeMult", "Mat");
2502     PetscCall(PetscOptionsEList("-matmattransmult_mpidense_mpidense_via", "Algorithmic approach", "MatMatTransposeMult", algTypes, nalg, algTypes[alg], &alg, &flg));
2503     PetscOptionsEnd();
2504   } else {
2505     PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatProduct_ABt", "Mat");
2506     PetscCall(PetscOptionsEList("-mat_product_algorithm", "Algorithmic approach", "MatProduct_ABt", algTypes, nalg, algTypes[alg], &alg, &flg));
2507     PetscOptionsEnd();
2508   }
2509   if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2510 
2511   C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_MPIDense_MPIDense;
2512   C->ops->productsymbolic          = MatProductSymbolic_ABt;
2513   PetscFunctionReturn(PETSC_SUCCESS);
2514 }
2515 
2516 static PetscErrorCode MatProductSetFromOptions_MPIDense(Mat C)
2517 {
2518   Mat_Product *product = C->product;
2519 
2520   PetscFunctionBegin;
2521   switch (product->type) {
2522   case MATPRODUCT_AB:
2523     PetscCall(MatProductSetFromOptions_MPIDense_AB(C));
2524     break;
2525   case MATPRODUCT_AtB:
2526     PetscCall(MatProductSetFromOptions_MPIDense_AtB(C));
2527     break;
2528   case MATPRODUCT_ABt:
2529     PetscCall(MatProductSetFromOptions_MPIDense_ABt(C));
2530     break;
2531   default:
2532     break;
2533   }
2534   PetscFunctionReturn(PETSC_SUCCESS);
2535 }
2536