xref: /petsc/src/mat/impls/is/matis.c (revision 174dc0c8cee294b82b85e4dd3b331b29396264fc)
1 /*
2     Creates a matrix class for using the Neumann-Neumann type preconditioners.
3     This stores the matrices in globally unassembled form. Each processor
4     assembles only its local Neumann problem and the parallel matrix vector
5     product is handled "implicitly".
6 
7     Currently this allows for only one subdomain per processor.
8 */
9 
10 #include <petsc/private/matisimpl.h> /*I "petscmat.h" I*/
11 #include <../src/mat/impls/aij/mpi/mpiaij.h>
12 #include <petsc/private/sfimpl.h>
13 #include <petsc/private/vecimpl.h>
14 #include <petsc/private/hashseti.h>
15 
16 #define MATIS_MAX_ENTRIES_INSERTION 2048
17 static PetscErrorCode MatSetValuesLocal_IS(Mat, PetscInt, const PetscInt *, PetscInt, const PetscInt *, const PetscScalar *, InsertMode);
18 static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat, PetscInt, const PetscInt *, PetscInt, const PetscInt *, const PetscScalar *, InsertMode);
19 static PetscErrorCode MatISSetUpScatters_Private(Mat);
20 
21 static PetscErrorCode MatISContainerDestroyPtAP_Private(void **ptr)
22 {
23   MatISPtAP ptap = (MatISPtAP)*ptr;
24 
25   PetscFunctionBegin;
26   PetscCall(MatDestroySubMatrices(ptap->ris1 ? 2 : 1, &ptap->lP));
27   PetscCall(ISDestroy(&ptap->cis0));
28   PetscCall(ISDestroy(&ptap->cis1));
29   PetscCall(ISDestroy(&ptap->ris0));
30   PetscCall(ISDestroy(&ptap->ris1));
31   PetscCall(PetscFree(ptap));
32   PetscFunctionReturn(PETSC_SUCCESS);
33 }
34 
35 static PetscErrorCode MatPtAPNumeric_IS_XAIJ(Mat A, Mat P, Mat C)
36 {
37   MatISPtAP      ptap;
38   Mat_IS        *matis = (Mat_IS *)A->data;
39   Mat            lA, lC;
40   MatReuse       reuse;
41   IS             ris[2], cis[2];
42   PetscContainer c;
43   PetscInt       n;
44 
45   PetscFunctionBegin;
46   PetscCall(PetscObjectQuery((PetscObject)C, "_MatIS_PtAP", (PetscObject *)&c));
47   PetscCheck(c, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing PtAP information");
48   PetscCall(PetscContainerGetPointer(c, (void **)&ptap));
49   ris[0] = ptap->ris0;
50   ris[1] = ptap->ris1;
51   cis[0] = ptap->cis0;
52   cis[1] = ptap->cis1;
53   n      = ptap->ris1 ? 2 : 1;
54   reuse  = ptap->lP ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
55   PetscCall(MatCreateSubMatrices(P, n, ris, cis, reuse, &ptap->lP));
56 
57   PetscCall(MatISGetLocalMat(A, &lA));
58   PetscCall(MatISGetLocalMat(C, &lC));
59   if (ptap->ris1) { /* unsymmetric A mapping */
60     Mat lPt;
61 
62     PetscCall(MatTranspose(ptap->lP[1], MAT_INITIAL_MATRIX, &lPt));
63     PetscCall(MatMatMatMult(lPt, lA, ptap->lP[0], reuse, ptap->fill, &lC));
64     if (matis->storel2l) PetscCall(PetscObjectCompose((PetscObject)A, "_MatIS_PtAP_l2l", (PetscObject)lPt));
65     PetscCall(MatDestroy(&lPt));
66   } else {
67     PetscCall(MatPtAP(lA, ptap->lP[0], reuse, ptap->fill, &lC));
68     if (matis->storel2l) PetscCall(PetscObjectCompose((PetscObject)C, "_MatIS_PtAP_l2l", (PetscObject)ptap->lP[0]));
69   }
70   if (reuse == MAT_INITIAL_MATRIX) {
71     PetscCall(MatISSetLocalMat(C, lC));
72     PetscCall(MatDestroy(&lC));
73   }
74   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
75   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
76   PetscFunctionReturn(PETSC_SUCCESS);
77 }
78 
79 static PetscErrorCode MatGetNonzeroColumnsLocal_Private(Mat PT, IS *cis)
80 {
81   Mat             Po, Pd;
82   IS              zd, zo;
83   const PetscInt *garray;
84   PetscInt       *aux, i, bs;
85   PetscInt        dc, stc, oc, ctd, cto;
86   PetscBool       ismpiaij, ismpibaij, isseqaij, isseqbaij;
87   MPI_Comm        comm;
88 
89   PetscFunctionBegin;
90   PetscValidHeaderSpecific(PT, MAT_CLASSID, 1);
91   PetscAssertPointer(cis, 2);
92   PetscCall(PetscObjectGetComm((PetscObject)PT, &comm));
93   bs = 1;
94   PetscCall(PetscObjectBaseTypeCompare((PetscObject)PT, MATMPIAIJ, &ismpiaij));
95   PetscCall(PetscObjectBaseTypeCompare((PetscObject)PT, MATMPIBAIJ, &ismpibaij));
96   PetscCall(PetscObjectBaseTypeCompare((PetscObject)PT, MATSEQAIJ, &isseqaij));
97   PetscCall(PetscObjectTypeCompare((PetscObject)PT, MATSEQBAIJ, &isseqbaij));
98   if (isseqaij || isseqbaij) {
99     Pd     = PT;
100     Po     = NULL;
101     garray = NULL;
102   } else if (ismpiaij) {
103     PetscCall(MatMPIAIJGetSeqAIJ(PT, &Pd, &Po, &garray));
104   } else if (ismpibaij) {
105     PetscCall(MatMPIBAIJGetSeqBAIJ(PT, &Pd, &Po, &garray));
106     PetscCall(MatGetBlockSize(PT, &bs));
107   } else SETERRQ(comm, PETSC_ERR_SUP, "Not for matrix type %s", ((PetscObject)PT)->type_name);
108 
109   /* identify any null columns in Pd or Po */
110   /* We use a tolerance comparison since it may happen that, with geometric multigrid,
111      some of the columns are not really zero, but very close to */
112   zo = zd = NULL;
113   if (Po) PetscCall(MatFindNonzeroRowsOrCols_Basic(Po, PETSC_TRUE, PETSC_SMALL, &zo));
114   PetscCall(MatFindNonzeroRowsOrCols_Basic(Pd, PETSC_TRUE, PETSC_SMALL, &zd));
115 
116   PetscCall(MatGetLocalSize(PT, NULL, &dc));
117   PetscCall(MatGetOwnershipRangeColumn(PT, &stc, NULL));
118   if (Po) PetscCall(MatGetLocalSize(Po, NULL, &oc));
119   else oc = 0;
120   PetscCall(PetscMalloc1((dc + oc) / bs, &aux));
121   if (zd) {
122     const PetscInt *idxs;
123     PetscInt        nz;
124 
125     /* this will throw an error if bs is not valid */
126     PetscCall(ISSetBlockSize(zd, bs));
127     PetscCall(ISGetLocalSize(zd, &nz));
128     PetscCall(ISGetIndices(zd, &idxs));
129     ctd = nz / bs;
130     for (i = 0; i < ctd; i++) aux[i] = (idxs[bs * i] + stc) / bs;
131     PetscCall(ISRestoreIndices(zd, &idxs));
132   } else {
133     ctd = dc / bs;
134     for (i = 0; i < ctd; i++) aux[i] = i + stc / bs;
135   }
136   if (zo) {
137     const PetscInt *idxs;
138     PetscInt        nz;
139 
140     /* this will throw an error if bs is not valid */
141     PetscCall(ISSetBlockSize(zo, bs));
142     PetscCall(ISGetLocalSize(zo, &nz));
143     PetscCall(ISGetIndices(zo, &idxs));
144     cto = nz / bs;
145     for (i = 0; i < cto; i++) aux[i + ctd] = garray[idxs[bs * i] / bs];
146     PetscCall(ISRestoreIndices(zo, &idxs));
147   } else {
148     cto = oc / bs;
149     for (i = 0; i < cto; i++) aux[i + ctd] = garray[i];
150   }
151   PetscCall(ISCreateBlock(comm, bs, ctd + cto, aux, PETSC_OWN_POINTER, cis));
152   PetscCall(ISDestroy(&zd));
153   PetscCall(ISDestroy(&zo));
154   PetscFunctionReturn(PETSC_SUCCESS);
155 }
156 
157 static PetscErrorCode MatPtAPSymbolic_IS_XAIJ(Mat A, Mat P, PetscReal fill, Mat C)
158 {
159   Mat                    PT, lA;
160   MatISPtAP              ptap;
161   ISLocalToGlobalMapping Crl2g, Ccl2g, rl2g, cl2g;
162   PetscContainer         c;
163   MatType                lmtype;
164   const PetscInt        *garray;
165   PetscInt               ibs, N, dc;
166   MPI_Comm               comm;
167 
168   PetscFunctionBegin;
169   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
170   PetscCall(MatSetType(C, MATIS));
171   PetscCall(MatISGetLocalMat(A, &lA));
172   PetscCall(MatGetType(lA, &lmtype));
173   PetscCall(MatISSetLocalMatType(C, lmtype));
174   PetscCall(MatGetSize(P, NULL, &N));
175   PetscCall(MatGetLocalSize(P, NULL, &dc));
176   PetscCall(MatSetSizes(C, dc, dc, N, N));
177   /* Not sure about this
178   PetscCall(MatGetBlockSizes(P,NULL,&ibs));
179   PetscCall(MatSetBlockSize(*C,ibs));
180 */
181 
182   PetscCall(PetscNew(&ptap));
183   PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &c));
184   PetscCall(PetscContainerSetPointer(c, ptap));
185   PetscCall(PetscContainerSetCtxDestroy(c, MatISContainerDestroyPtAP_Private));
186   PetscCall(PetscObjectCompose((PetscObject)C, "_MatIS_PtAP", (PetscObject)c));
187   PetscCall(PetscContainerDestroy(&c));
188   ptap->fill = fill;
189 
190   PetscCall(MatISGetLocalToGlobalMapping(A, &rl2g, &cl2g));
191 
192   PetscCall(ISLocalToGlobalMappingGetBlockSize(cl2g, &ibs));
193   PetscCall(ISLocalToGlobalMappingGetSize(cl2g, &N));
194   PetscCall(ISLocalToGlobalMappingGetBlockIndices(cl2g, &garray));
195   PetscCall(ISCreateBlock(comm, ibs, N / ibs, garray, PETSC_COPY_VALUES, &ptap->ris0));
196   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(cl2g, &garray));
197 
198   PetscCall(MatCreateSubMatrix(P, ptap->ris0, NULL, MAT_INITIAL_MATRIX, &PT));
199   PetscCall(MatGetNonzeroColumnsLocal_Private(PT, &ptap->cis0));
200   PetscCall(ISLocalToGlobalMappingCreateIS(ptap->cis0, &Ccl2g));
201   PetscCall(MatDestroy(&PT));
202 
203   Crl2g = NULL;
204   if (rl2g != cl2g) { /* unsymmetric A mapping */
205     PetscBool same, lsame = PETSC_FALSE;
206     PetscInt  N1, ibs1;
207 
208     PetscCall(ISLocalToGlobalMappingGetSize(rl2g, &N1));
209     PetscCall(ISLocalToGlobalMappingGetBlockSize(rl2g, &ibs1));
210     PetscCall(ISLocalToGlobalMappingGetBlockIndices(rl2g, &garray));
211     PetscCall(ISCreateBlock(comm, ibs, N1 / ibs, garray, PETSC_COPY_VALUES, &ptap->ris1));
212     PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(rl2g, &garray));
213     if (ibs1 == ibs && N1 == N) { /* check if the l2gmaps are the same */
214       const PetscInt *i1, *i2;
215 
216       PetscCall(ISBlockGetIndices(ptap->ris0, &i1));
217       PetscCall(ISBlockGetIndices(ptap->ris1, &i2));
218       PetscCall(PetscArraycmp(i1, i2, N, &lsame));
219     }
220     PetscCallMPI(MPIU_Allreduce(&lsame, &same, 1, MPIU_BOOL, MPI_LAND, comm));
221     if (same) {
222       PetscCall(ISDestroy(&ptap->ris1));
223     } else {
224       PetscCall(MatCreateSubMatrix(P, ptap->ris1, NULL, MAT_INITIAL_MATRIX, &PT));
225       PetscCall(MatGetNonzeroColumnsLocal_Private(PT, &ptap->cis1));
226       PetscCall(ISLocalToGlobalMappingCreateIS(ptap->cis1, &Crl2g));
227       PetscCall(MatDestroy(&PT));
228     }
229   }
230   /* Not sure about this
231   if (!Crl2g) {
232     PetscCall(MatGetBlockSize(C,&ibs));
233     PetscCall(ISLocalToGlobalMappingSetBlockSize(Ccl2g,ibs));
234   }
235 */
236   PetscCall(MatSetLocalToGlobalMapping(C, Crl2g ? Crl2g : Ccl2g, Ccl2g));
237   PetscCall(ISLocalToGlobalMappingDestroy(&Crl2g));
238   PetscCall(ISLocalToGlobalMappingDestroy(&Ccl2g));
239 
240   C->ops->ptapnumeric = MatPtAPNumeric_IS_XAIJ;
241   PetscFunctionReturn(PETSC_SUCCESS);
242 }
243 
244 static PetscErrorCode MatProductSymbolic_PtAP_IS_XAIJ(Mat C)
245 {
246   Mat_Product *product = C->product;
247   Mat          A = product->A, P = product->B;
248   PetscReal    fill = product->fill;
249 
250   PetscFunctionBegin;
251   PetscCall(MatPtAPSymbolic_IS_XAIJ(A, P, fill, C));
252   C->ops->productnumeric = MatProductNumeric_PtAP;
253   PetscFunctionReturn(PETSC_SUCCESS);
254 }
255 
256 static PetscErrorCode MatProductSetFromOptions_IS_XAIJ_PtAP(Mat C)
257 {
258   PetscFunctionBegin;
259   C->ops->productsymbolic = MatProductSymbolic_PtAP_IS_XAIJ;
260   PetscFunctionReturn(PETSC_SUCCESS);
261 }
262 
263 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_IS_XAIJ(Mat C)
264 {
265   Mat_Product *product = C->product;
266 
267   PetscFunctionBegin;
268   if (product->type == MATPRODUCT_PtAP) PetscCall(MatProductSetFromOptions_IS_XAIJ_PtAP(C));
269   PetscFunctionReturn(PETSC_SUCCESS);
270 }
271 
272 static PetscErrorCode MatISContainerDestroyFields_Private(void **ptr)
273 {
274   MatISLocalFields lf = (MatISLocalFields)*ptr;
275   PetscInt         i;
276 
277   PetscFunctionBegin;
278   for (i = 0; i < lf->nr; i++) PetscCall(ISDestroy(&lf->rf[i]));
279   for (i = 0; i < lf->nc; i++) PetscCall(ISDestroy(&lf->cf[i]));
280   PetscCall(PetscFree2(lf->rf, lf->cf));
281   PetscCall(PetscFree(lf));
282   PetscFunctionReturn(PETSC_SUCCESS);
283 }
284 
285 static PetscErrorCode MatConvert_SeqXAIJ_IS(Mat A, MatType type, MatReuse reuse, Mat *newmat)
286 {
287   Mat B, lB;
288 
289   PetscFunctionBegin;
290   if (reuse != MAT_REUSE_MATRIX) {
291     ISLocalToGlobalMapping rl2g, cl2g;
292     PetscInt               bs;
293     IS                     is;
294 
295     PetscCall(MatGetBlockSize(A, &bs));
296     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A), A->rmap->n / bs, 0, 1, &is));
297     if (bs > 1) {
298       IS       is2;
299       PetscInt i, *aux;
300 
301       PetscCall(ISGetLocalSize(is, &i));
302       PetscCall(ISGetIndices(is, (const PetscInt **)&aux));
303       PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A), bs, i, aux, PETSC_COPY_VALUES, &is2));
304       PetscCall(ISRestoreIndices(is, (const PetscInt **)&aux));
305       PetscCall(ISDestroy(&is));
306       is = is2;
307     }
308     PetscCall(ISSetIdentity(is));
309     PetscCall(ISLocalToGlobalMappingCreateIS(is, &rl2g));
310     PetscCall(ISDestroy(&is));
311     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A), A->cmap->n / bs, 0, 1, &is));
312     if (bs > 1) {
313       IS       is2;
314       PetscInt i, *aux;
315 
316       PetscCall(ISGetLocalSize(is, &i));
317       PetscCall(ISGetIndices(is, (const PetscInt **)&aux));
318       PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A), bs, i, aux, PETSC_COPY_VALUES, &is2));
319       PetscCall(ISRestoreIndices(is, (const PetscInt **)&aux));
320       PetscCall(ISDestroy(&is));
321       is = is2;
322     }
323     PetscCall(ISSetIdentity(is));
324     PetscCall(ISLocalToGlobalMappingCreateIS(is, &cl2g));
325     PetscCall(ISDestroy(&is));
326     PetscCall(MatCreateIS(PetscObjectComm((PetscObject)A), bs, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N, rl2g, cl2g, &B));
327     PetscCall(ISLocalToGlobalMappingDestroy(&rl2g));
328     PetscCall(ISLocalToGlobalMappingDestroy(&cl2g));
329     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &lB));
330     if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
331   } else {
332     B = *newmat;
333     PetscCall(PetscObjectReference((PetscObject)A));
334     lB = A;
335   }
336   PetscCall(MatISSetLocalMat(B, lB));
337   PetscCall(MatDestroy(&lB));
338   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
339   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
340   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
341   PetscFunctionReturn(PETSC_SUCCESS);
342 }
343 
344 static PetscErrorCode MatISScaleDisassembling_Private(Mat A)
345 {
346   Mat_IS         *matis = (Mat_IS *)A->data;
347   PetscScalar    *aa;
348   const PetscInt *ii, *jj;
349   PetscInt        i, n, m;
350   PetscInt       *ecount, **eneighs;
351   PetscBool       flg;
352 
353   PetscFunctionBegin;
354   PetscCall(MatGetRowIJ(matis->A, 0, PETSC_FALSE, PETSC_FALSE, &m, &ii, &jj, &flg));
355   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot get IJ structure");
356   PetscCall(ISLocalToGlobalMappingGetNodeInfo(matis->rmapping, &n, &ecount, &eneighs));
357   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected %" PetscInt_FMT " != %" PetscInt_FMT, m, n);
358   PetscCall(MatSeqAIJGetArray(matis->A, &aa));
359   for (i = 0; i < n; i++) {
360     if (ecount[i] > 1) {
361       PetscInt j;
362 
363       for (j = ii[i]; j < ii[i + 1]; j++) {
364         PetscInt  i2   = jj[j], p, p2;
365         PetscReal scal = 0.0;
366 
367         for (p = 0; p < ecount[i]; p++) {
368           for (p2 = 0; p2 < ecount[i2]; p2++) {
369             if (eneighs[i][p] == eneighs[i2][p2]) {
370               scal += 1.0;
371               break;
372             }
373           }
374         }
375         if (scal) aa[j] /= scal;
376       }
377     }
378   }
379   PetscCall(ISLocalToGlobalMappingRestoreNodeInfo(matis->rmapping, &n, &ecount, &eneighs));
380   PetscCall(MatSeqAIJRestoreArray(matis->A, &aa));
381   PetscCall(MatRestoreRowIJ(matis->A, 0, PETSC_FALSE, PETSC_FALSE, &m, &ii, &jj, &flg));
382   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot restore IJ structure");
383   PetscFunctionReturn(PETSC_SUCCESS);
384 }
385 
386 typedef enum {
387   MAT_IS_DISASSEMBLE_L2G_NATURAL,
388   MAT_IS_DISASSEMBLE_L2G_MAT,
389   MAT_IS_DISASSEMBLE_L2G_ND
390 } MatISDisassemblel2gType;
391 
392 static PetscErrorCode MatMPIXAIJComputeLocalToGlobalMapping_Private(Mat A, ISLocalToGlobalMapping *l2g)
393 {
394   Mat                     Ad, Ao;
395   IS                      is, ndmap, ndsub;
396   MPI_Comm                comm;
397   const PetscInt         *garray, *ndmapi;
398   PetscInt                bs, i, cnt, nl, *ncount, *ndmapc;
399   PetscBool               ismpiaij, ismpibaij;
400   const char *const       MatISDisassemblel2gTypes[] = {"NATURAL", "MAT", "ND", "MatISDisassemblel2gType", "MAT_IS_DISASSEMBLE_L2G_", NULL};
401   MatISDisassemblel2gType mode                       = MAT_IS_DISASSEMBLE_L2G_NATURAL;
402   MatPartitioning         part;
403   PetscSF                 sf;
404   PetscObject             dm;
405 
406   PetscFunctionBegin;
407   PetscOptionsBegin(PetscObjectComm((PetscObject)A), ((PetscObject)A)->prefix, "MatIS l2g disassembling options", "Mat");
408   PetscCall(PetscOptionsEnum("-mat_is_disassemble_l2g_type", "Type of local-to-global mapping to be used for disassembling", "MatISDisassemblel2gType", MatISDisassemblel2gTypes, (PetscEnum)mode, (PetscEnum *)&mode, NULL));
409   PetscOptionsEnd();
410   if (mode == MAT_IS_DISASSEMBLE_L2G_MAT) {
411     PetscCall(MatGetLocalToGlobalMapping(A, l2g, NULL));
412     PetscFunctionReturn(PETSC_SUCCESS);
413   }
414   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
415   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &ismpiaij));
416   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIBAIJ, &ismpibaij));
417   PetscCall(MatGetBlockSize(A, &bs));
418   switch (mode) {
419   case MAT_IS_DISASSEMBLE_L2G_ND:
420     PetscCall(MatPartitioningCreate(comm, &part));
421     PetscCall(MatPartitioningSetAdjacency(part, A));
422     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part, ((PetscObject)A)->prefix));
423     PetscCall(MatPartitioningSetFromOptions(part));
424     PetscCall(MatPartitioningApplyND(part, &ndmap));
425     PetscCall(MatPartitioningDestroy(&part));
426     PetscCall(ISBuildTwoSided(ndmap, NULL, &ndsub));
427     PetscCall(MatMPIAIJSetUseScalableIncreaseOverlap(A, PETSC_TRUE));
428     PetscCall(MatIncreaseOverlap(A, 1, &ndsub, 1));
429     PetscCall(ISLocalToGlobalMappingCreateIS(ndsub, l2g));
430 
431     /* it may happen that a separator node is not properly shared */
432     PetscCall(ISLocalToGlobalMappingGetNodeInfo(*l2g, &nl, &ncount, NULL));
433     PetscCall(PetscSFCreate(comm, &sf));
434     PetscCall(ISLocalToGlobalMappingGetIndices(*l2g, &garray));
435     PetscCall(PetscSFSetGraphLayout(sf, A->rmap, nl, NULL, PETSC_OWN_POINTER, garray));
436     PetscCall(ISLocalToGlobalMappingRestoreIndices(*l2g, &garray));
437     PetscCall(PetscCalloc1(A->rmap->n, &ndmapc));
438     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, ncount, ndmapc, MPI_REPLACE));
439     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, ncount, ndmapc, MPI_REPLACE));
440     PetscCall(ISLocalToGlobalMappingRestoreNodeInfo(*l2g, NULL, &ncount, NULL));
441     PetscCall(ISGetIndices(ndmap, &ndmapi));
442     for (i = 0, cnt = 0; i < A->rmap->n; i++)
443       if (ndmapi[i] < 0 && ndmapc[i] < 2) cnt++;
444 
445     PetscCallMPI(MPIU_Allreduce(&cnt, &i, 1, MPIU_INT, MPI_MAX, comm));
446     if (i) { /* we detected isolated separator nodes */
447       Mat                    A2, A3;
448       IS                    *workis, is2;
449       PetscScalar           *vals;
450       PetscInt               gcnt = i, *dnz, *onz, j, *lndmapi;
451       ISLocalToGlobalMapping ll2g;
452       PetscBool              flg;
453       const PetscInt        *ii, *jj;
454 
455       /* communicate global id of separators */
456       MatPreallocateBegin(comm, A->rmap->n, A->cmap->n, dnz, onz);
457       for (i = 0, cnt = 0; i < A->rmap->n; i++) dnz[i] = ndmapi[i] < 0 ? i + A->rmap->rstart : -1;
458 
459       PetscCall(PetscMalloc1(nl, &lndmapi));
460       PetscCall(PetscSFBcastBegin(sf, MPIU_INT, dnz, lndmapi, MPI_REPLACE));
461 
462       /* compute adjacency of isolated separators node */
463       PetscCall(PetscMalloc1(gcnt, &workis));
464       for (i = 0, cnt = 0; i < A->rmap->n; i++) {
465         if (ndmapi[i] < 0 && ndmapc[i] < 2) PetscCall(ISCreateStride(comm, 1, i + A->rmap->rstart, 1, &workis[cnt++]));
466       }
467       for (i = cnt; i < gcnt; i++) PetscCall(ISCreateStride(comm, 0, 0, 1, &workis[i]));
468       for (i = 0; i < gcnt; i++) {
469         PetscCall(PetscObjectSetName((PetscObject)workis[i], "ISOLATED"));
470         PetscCall(ISViewFromOptions(workis[i], NULL, "-view_isolated_separators"));
471       }
472 
473       /* no communications since all the ISes correspond to locally owned rows */
474       PetscCall(MatIncreaseOverlap(A, gcnt, workis, 1));
475 
476       /* end communicate global id of separators */
477       PetscCall(PetscSFBcastEnd(sf, MPIU_INT, dnz, lndmapi, MPI_REPLACE));
478 
479       /* communicate new layers : create a matrix and transpose it */
480       PetscCall(PetscArrayzero(dnz, A->rmap->n));
481       PetscCall(PetscArrayzero(onz, A->rmap->n));
482       for (i = 0, j = 0; i < A->rmap->n; i++) {
483         if (ndmapi[i] < 0 && ndmapc[i] < 2) {
484           const PetscInt *idxs;
485           PetscInt        s;
486 
487           PetscCall(ISGetLocalSize(workis[j], &s));
488           PetscCall(ISGetIndices(workis[j], &idxs));
489           PetscCall(MatPreallocateSet(i + A->rmap->rstart, s, idxs, dnz, onz));
490           j++;
491         }
492       }
493       PetscCheck(j == cnt, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected local count %" PetscInt_FMT " != %" PetscInt_FMT, j, cnt);
494 
495       for (i = 0; i < gcnt; i++) {
496         PetscCall(PetscObjectSetName((PetscObject)workis[i], "EXTENDED"));
497         PetscCall(ISViewFromOptions(workis[i], NULL, "-view_isolated_separators"));
498       }
499 
500       for (i = 0, j = 0; i < A->rmap->n; i++) j = PetscMax(j, dnz[i] + onz[i]);
501       PetscCall(PetscMalloc1(j, &vals));
502       for (i = 0; i < j; i++) vals[i] = 1.0;
503 
504       PetscCall(MatCreate(comm, &A2));
505       PetscCall(MatSetType(A2, MATMPIAIJ));
506       PetscCall(MatSetSizes(A2, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
507       PetscCall(MatMPIAIJSetPreallocation(A2, 0, dnz, 0, onz));
508       PetscCall(MatSetOption(A2, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
509       for (i = 0, j = 0; i < A2->rmap->n; i++) {
510         PetscInt        row = i + A2->rmap->rstart, s = dnz[i] + onz[i];
511         const PetscInt *idxs;
512 
513         if (s) {
514           PetscCall(ISGetIndices(workis[j], &idxs));
515           PetscCall(MatSetValues(A2, 1, &row, s, idxs, vals, INSERT_VALUES));
516           PetscCall(ISRestoreIndices(workis[j], &idxs));
517           j++;
518         }
519       }
520       PetscCheck(j == cnt, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected local count %" PetscInt_FMT " != %" PetscInt_FMT, j, cnt);
521       PetscCall(PetscFree(vals));
522       PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY));
523       PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY));
524       PetscCall(MatTranspose(A2, MAT_INPLACE_MATRIX, &A2));
525 
526       /* extract submatrix corresponding to the coupling "owned separators" x "isolated separators" */
527       for (i = 0, j = 0; i < nl; i++)
528         if (lndmapi[i] >= 0) lndmapi[j++] = lndmapi[i];
529       PetscCall(ISCreateGeneral(comm, j, lndmapi, PETSC_USE_POINTER, &is));
530       PetscCall(MatMPIAIJGetLocalMatCondensed(A2, MAT_INITIAL_MATRIX, &is, NULL, &A3));
531       PetscCall(ISDestroy(&is));
532       PetscCall(MatDestroy(&A2));
533 
534       /* extend local to global map to include connected isolated separators */
535       PetscCall(PetscObjectQuery((PetscObject)A3, "_petsc_GetLocalMatCondensed_iscol", (PetscObject *)&is));
536       PetscCheck(is, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing column map");
537       PetscCall(ISLocalToGlobalMappingCreateIS(is, &ll2g));
538       PetscCall(MatGetRowIJ(A3, 0, PETSC_FALSE, PETSC_FALSE, &i, &ii, &jj, &flg));
539       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot get IJ structure");
540       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ii[i], jj, PETSC_COPY_VALUES, &is));
541       PetscCall(MatRestoreRowIJ(A3, 0, PETSC_FALSE, PETSC_FALSE, &i, &ii, &jj, &flg));
542       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot get IJ structure");
543       PetscCall(ISLocalToGlobalMappingApplyIS(ll2g, is, &is2));
544       PetscCall(ISDestroy(&is));
545       PetscCall(ISLocalToGlobalMappingDestroy(&ll2g));
546 
547       /* add new nodes to the local-to-global map */
548       PetscCall(ISLocalToGlobalMappingDestroy(l2g));
549       PetscCall(ISExpand(ndsub, is2, &is));
550       PetscCall(ISDestroy(&is2));
551       PetscCall(ISLocalToGlobalMappingCreateIS(is, l2g));
552       PetscCall(ISDestroy(&is));
553 
554       PetscCall(MatDestroy(&A3));
555       PetscCall(PetscFree(lndmapi));
556       MatPreallocateEnd(dnz, onz);
557       for (i = 0; i < gcnt; i++) PetscCall(ISDestroy(&workis[i]));
558       PetscCall(PetscFree(workis));
559     }
560     PetscCall(ISRestoreIndices(ndmap, &ndmapi));
561     PetscCall(PetscSFDestroy(&sf));
562     PetscCall(PetscFree(ndmapc));
563     PetscCall(ISDestroy(&ndmap));
564     PetscCall(ISDestroy(&ndsub));
565     PetscCall(ISLocalToGlobalMappingSetBlockSize(*l2g, bs));
566     PetscCall(ISLocalToGlobalMappingViewFromOptions(*l2g, NULL, "-mat_is_nd_l2g_view"));
567     break;
568   case MAT_IS_DISASSEMBLE_L2G_NATURAL:
569     PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_dm", &dm));
570     if (dm) { /* if a matrix comes from a DM, most likely we can use the l2gmap if any */
571       PetscCall(MatGetLocalToGlobalMapping(A, l2g, NULL));
572       PetscCall(PetscObjectReference((PetscObject)*l2g));
573       if (*l2g) PetscFunctionReturn(PETSC_SUCCESS);
574     }
575     if (ismpiaij) {
576       PetscCall(MatMPIAIJGetSeqAIJ(A, &Ad, &Ao, &garray));
577     } else if (ismpibaij) {
578       PetscCall(MatMPIBAIJGetSeqBAIJ(A, &Ad, &Ao, &garray));
579     } else SETERRQ(comm, PETSC_ERR_SUP, "Type %s", ((PetscObject)A)->type_name);
580     if (A->rmap->n) {
581       PetscInt dc, oc, stc, *aux;
582 
583       PetscCall(MatGetLocalSize(Ad, NULL, &dc));
584       PetscCall(MatGetLocalSize(Ao, NULL, &oc));
585       PetscCheck(!oc || garray, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "garray not present");
586       PetscCall(MatGetOwnershipRangeColumn(A, &stc, NULL));
587       PetscCall(PetscMalloc1((dc + oc) / bs, &aux));
588       for (i = 0; i < dc / bs; i++) aux[i] = i + stc / bs;
589       for (i = 0; i < oc / bs; i++) aux[i + dc / bs] = (ismpiaij ? garray[i * bs] / bs : garray[i]);
590       PetscCall(ISCreateBlock(comm, bs, (dc + oc) / bs, aux, PETSC_OWN_POINTER, &is));
591     } else {
592       PetscCall(ISCreateBlock(comm, 1, 0, NULL, PETSC_OWN_POINTER, &is));
593     }
594     PetscCall(ISLocalToGlobalMappingCreateIS(is, l2g));
595     PetscCall(ISDestroy(&is));
596     break;
597   default:
598     SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unsupported l2g disassembling type %d", mode);
599   }
600   PetscFunctionReturn(PETSC_SUCCESS);
601 }
602 
603 PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat A, MatType type, MatReuse reuse, Mat *newmat)
604 {
605   Mat                    lA, Ad, Ao, B = NULL;
606   ISLocalToGlobalMapping rl2g, cl2g;
607   IS                     is;
608   MPI_Comm               comm;
609   void                  *ptrs[2];
610   const char            *names[2] = {"_convert_csr_aux", "_convert_csr_data"};
611   const PetscInt        *garray;
612   PetscScalar           *dd, *od, *aa, *data;
613   const PetscInt        *di, *dj, *oi, *oj;
614   const PetscInt        *odi, *odj, *ooi, *ooj;
615   PetscInt              *aux, *ii, *jj;
616   PetscInt               rbs, cbs, lc, dr, dc, oc, str, stc, nnz, i, jd, jo, cum;
617   PetscBool              flg, ismpiaij, ismpibaij, was_inplace = PETSC_FALSE, cong;
618   PetscMPIInt            size;
619 
620   PetscFunctionBegin;
621   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
622   PetscCallMPI(MPI_Comm_size(comm, &size));
623   if (size == 1) {
624     PetscCall(MatConvert_SeqXAIJ_IS(A, type, reuse, newmat));
625     PetscFunctionReturn(PETSC_SUCCESS);
626   }
627   PetscCall(MatGetBlockSizes(A, &rbs, &cbs));
628   PetscCall(MatHasCongruentLayouts(A, &cong));
629   if (reuse != MAT_REUSE_MATRIX && cong && rbs == cbs) {
630     PetscCall(MatMPIXAIJComputeLocalToGlobalMapping_Private(A, &rl2g));
631     PetscCall(MatCreate(comm, &B));
632     PetscCall(MatSetType(B, MATIS));
633     PetscCall(MatSetSizes(B, A->rmap->n, A->rmap->n, A->rmap->N, A->rmap->N));
634     PetscCall(MatSetLocalToGlobalMapping(B, rl2g, rl2g));
635     PetscCall(MatSetBlockSizes(B, rbs, rbs));
636     PetscCall(ISLocalToGlobalMappingDestroy(&rl2g));
637     if (reuse == MAT_INPLACE_MATRIX) was_inplace = PETSC_TRUE;
638     reuse = MAT_REUSE_MATRIX;
639   }
640   if (reuse == MAT_REUSE_MATRIX) {
641     Mat            *newlA, lA;
642     IS              rows, cols;
643     const PetscInt *ridx, *cidx;
644     PetscInt        nr, nc;
645 
646     if (!B) B = *newmat;
647     PetscCall(MatISGetLocalToGlobalMapping(B, &rl2g, &cl2g));
648     PetscCall(ISLocalToGlobalMappingGetBlockIndices(rl2g, &ridx));
649     PetscCall(ISLocalToGlobalMappingGetBlockIndices(cl2g, &cidx));
650     PetscCall(ISLocalToGlobalMappingGetSize(rl2g, &nr));
651     PetscCall(ISLocalToGlobalMappingGetSize(cl2g, &nc));
652     PetscCall(ISLocalToGlobalMappingGetBlockSize(rl2g, &rbs));
653     PetscCall(ISLocalToGlobalMappingGetBlockSize(cl2g, &cbs));
654     PetscCall(ISCreateBlock(comm, rbs, nr / rbs, ridx, PETSC_USE_POINTER, &rows));
655     if (rl2g != cl2g) {
656       PetscCall(ISCreateBlock(comm, cbs, nc / cbs, cidx, PETSC_USE_POINTER, &cols));
657     } else {
658       PetscCall(PetscObjectReference((PetscObject)rows));
659       cols = rows;
660     }
661     PetscCall(MatISGetLocalMat(B, &lA));
662     PetscCall(MatCreateSubMatrices(A, 1, &rows, &cols, MAT_INITIAL_MATRIX, &newlA));
663     PetscCall(MatConvert(newlA[0], MATSEQAIJ, MAT_INPLACE_MATRIX, &newlA[0]));
664     PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(rl2g, &ridx));
665     PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(cl2g, &cidx));
666     PetscCall(ISDestroy(&rows));
667     PetscCall(ISDestroy(&cols));
668     if (!lA->preallocated) { /* first time */
669       PetscCall(MatDuplicate(newlA[0], MAT_COPY_VALUES, &lA));
670       PetscCall(MatISSetLocalMat(B, lA));
671       PetscCall(PetscObjectDereference((PetscObject)lA));
672     }
673     PetscCall(MatCopy(newlA[0], lA, SAME_NONZERO_PATTERN));
674     PetscCall(MatDestroySubMatrices(1, &newlA));
675     PetscCall(MatISScaleDisassembling_Private(B));
676     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
677     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
678     if (was_inplace) PetscCall(MatHeaderReplace(A, &B));
679     else *newmat = B;
680     PetscFunctionReturn(PETSC_SUCCESS);
681   }
682   /* general case, just compress out the column space */
683   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &ismpiaij));
684   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIBAIJ, &ismpibaij));
685   if (ismpiaij) {
686     cbs = 1; /* We cannot guarantee the off-process matrix will respect the column block size */
687     PetscCall(MatMPIAIJGetSeqAIJ(A, &Ad, &Ao, &garray));
688   } else if (ismpibaij) {
689     PetscCall(MatMPIBAIJGetSeqBAIJ(A, &Ad, &Ao, &garray));
690     PetscCall(MatConvert(Ad, MATSEQAIJ, MAT_INITIAL_MATRIX, &Ad));
691     PetscCall(MatConvert(Ao, MATSEQAIJ, MAT_INITIAL_MATRIX, &Ao));
692   } else SETERRQ(comm, PETSC_ERR_SUP, "Type %s", ((PetscObject)A)->type_name);
693   PetscCall(MatSeqAIJGetArray(Ad, &dd));
694   PetscCall(MatSeqAIJGetArray(Ao, &od));
695 
696   /* access relevant information from MPIAIJ */
697   PetscCall(MatGetOwnershipRange(A, &str, NULL));
698   PetscCall(MatGetOwnershipRangeColumn(A, &stc, NULL));
699   PetscCall(MatGetLocalSize(Ad, &dr, &dc));
700   PetscCall(MatGetLocalSize(Ao, NULL, &oc));
701   PetscCheck(!oc || garray, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "garray not present");
702 
703   PetscCall(MatGetRowIJ(Ad, 0, PETSC_FALSE, PETSC_FALSE, &i, &di, &dj, &flg));
704   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot get IJ structure");
705   PetscCall(MatGetRowIJ(Ao, 0, PETSC_FALSE, PETSC_FALSE, &i, &oi, &oj, &flg));
706   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot get IJ structure");
707   nnz = di[dr] + oi[dr];
708   /* store original pointers to be restored later */
709   odi = di;
710   odj = dj;
711   ooi = oi;
712   ooj = oj;
713 
714   /* generate l2g maps for rows and cols */
715   PetscCall(ISCreateStride(comm, dr / rbs, str / rbs, 1, &is));
716   if (rbs > 1) {
717     IS is2;
718 
719     PetscCall(ISGetLocalSize(is, &i));
720     PetscCall(ISGetIndices(is, (const PetscInt **)&aux));
721     PetscCall(ISCreateBlock(comm, rbs, i, aux, PETSC_COPY_VALUES, &is2));
722     PetscCall(ISRestoreIndices(is, (const PetscInt **)&aux));
723     PetscCall(ISDestroy(&is));
724     is = is2;
725   }
726   PetscCall(ISLocalToGlobalMappingCreateIS(is, &rl2g));
727   PetscCall(ISDestroy(&is));
728   if (dr) {
729     PetscCall(PetscMalloc1((dc + oc) / cbs, &aux));
730     for (i = 0; i < dc / cbs; i++) aux[i] = i + stc / cbs;
731     for (i = 0; i < oc / cbs; i++) aux[i + dc / cbs] = garray[i];
732     PetscCall(ISCreateBlock(comm, cbs, (dc + oc) / cbs, aux, PETSC_OWN_POINTER, &is));
733     lc = dc + oc;
734   } else {
735     PetscCall(ISCreateBlock(comm, cbs, 0, NULL, PETSC_OWN_POINTER, &is));
736     lc = 0;
737   }
738   PetscCall(ISLocalToGlobalMappingCreateIS(is, &cl2g));
739   PetscCall(ISDestroy(&is));
740 
741   /* create MATIS object */
742   PetscCall(MatCreate(comm, &B));
743   PetscCall(MatSetSizes(B, dr, dc, PETSC_DECIDE, PETSC_DECIDE));
744   PetscCall(MatSetType(B, MATIS));
745   PetscCall(MatSetBlockSizes(B, rbs, cbs));
746   PetscCall(MatSetLocalToGlobalMapping(B, rl2g, cl2g));
747   PetscCall(ISLocalToGlobalMappingDestroy(&rl2g));
748   PetscCall(ISLocalToGlobalMappingDestroy(&cl2g));
749 
750   /* merge local matrices */
751   PetscCall(PetscMalloc1(nnz + dr + 1, &aux));
752   PetscCall(PetscMalloc1(nnz, &data));
753   ii  = aux;
754   jj  = aux + dr + 1;
755   aa  = data;
756   *ii = *(di++) + *(oi++);
757   for (jd = 0, jo = 0, cum = 0; *ii < nnz; cum++) {
758     for (; jd < *di; jd++) {
759       *jj++ = *dj++;
760       *aa++ = *dd++;
761     }
762     for (; jo < *oi; jo++) {
763       *jj++ = *oj++ + dc;
764       *aa++ = *od++;
765     }
766     *(++ii) = *(di++) + *(oi++);
767   }
768   for (; cum < dr; cum++) *(++ii) = nnz;
769 
770   PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_FALSE, PETSC_FALSE, &i, &odi, &odj, &flg));
771   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot restore IJ structure");
772   PetscCall(MatRestoreRowIJ(Ao, 0, PETSC_FALSE, PETSC_FALSE, &i, &ooi, &ooj, &flg));
773   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot restore IJ structure");
774   PetscCall(MatSeqAIJRestoreArray(Ad, &dd));
775   PetscCall(MatSeqAIJRestoreArray(Ao, &od));
776 
777   ii = aux;
778   jj = aux + dr + 1;
779   aa = data;
780   PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, dr, lc, ii, jj, aa, &lA));
781 
782   /* create containers to destroy the data */
783   ptrs[0] = aux;
784   ptrs[1] = data;
785   for (i = 0; i < 2; i++) PetscCall(PetscObjectContainerCompose((PetscObject)lA, names[i], ptrs[i], PetscCtxDestroyDefault));
786   if (ismpibaij) { /* destroy converted local matrices */
787     PetscCall(MatDestroy(&Ad));
788     PetscCall(MatDestroy(&Ao));
789   }
790 
791   /* finalize matrix */
792   PetscCall(MatISSetLocalMat(B, lA));
793   PetscCall(MatDestroy(&lA));
794   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
795   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
796   if (reuse == MAT_INPLACE_MATRIX) {
797     PetscCall(MatHeaderReplace(A, &B));
798   } else *newmat = B;
799   PetscFunctionReturn(PETSC_SUCCESS);
800 }
801 
802 PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat A, MatType type, MatReuse reuse, Mat *newmat)
803 {
804   Mat                  **nest, *snest, **rnest, lA, B;
805   IS                    *iscol, *isrow, *islrow, *islcol;
806   ISLocalToGlobalMapping rl2g, cl2g;
807   MPI_Comm               comm;
808   PetscInt              *lr, *lc, *l2gidxs;
809   PetscInt               i, j, nr, nc, rbs, cbs;
810   PetscBool              convert, lreuse, *istrans;
811   PetscBool3             allow_repeated = PETSC_BOOL3_UNKNOWN;
812 
813   PetscFunctionBegin;
814   PetscCall(MatNestGetSubMats(A, &nr, &nc, &nest));
815   lreuse = PETSC_FALSE;
816   rnest  = NULL;
817   if (reuse == MAT_REUSE_MATRIX) {
818     PetscBool ismatis, isnest;
819 
820     PetscCall(PetscObjectTypeCompare((PetscObject)*newmat, MATIS, &ismatis));
821     PetscCheck(ismatis, PetscObjectComm((PetscObject)*newmat), PETSC_ERR_USER, "Cannot reuse matrix of type %s", ((PetscObject)*newmat)->type_name);
822     PetscCall(MatISGetLocalMat(*newmat, &lA));
823     PetscCall(PetscObjectTypeCompare((PetscObject)lA, MATNEST, &isnest));
824     if (isnest) {
825       PetscCall(MatNestGetSubMats(lA, &i, &j, &rnest));
826       lreuse = (PetscBool)(i == nr && j == nc);
827       if (!lreuse) rnest = NULL;
828     }
829   }
830   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
831   PetscCall(PetscCalloc2(nr, &lr, nc, &lc));
832   PetscCall(PetscCalloc6(nr, &isrow, nc, &iscol, nr, &islrow, nc, &islcol, nr * nc, &snest, nr * nc, &istrans));
833   PetscCall(MatNestGetISs(A, isrow, iscol));
834   for (i = 0; i < nr; i++) {
835     for (j = 0; j < nc; j++) {
836       PetscBool ismatis, sallow;
837       PetscInt  l1, l2, lb1, lb2, ij = i * nc + j;
838 
839       /* Null matrix pointers are allowed in MATNEST */
840       if (!nest[i][j]) continue;
841 
842       /* Nested matrices should be of type MATIS */
843       PetscCall(PetscObjectTypeCompare((PetscObject)nest[i][j], MATTRANSPOSEVIRTUAL, &istrans[ij]));
844       if (istrans[ij]) {
845         Mat T, lT;
846         PetscCall(MatTransposeGetMat(nest[i][j], &T));
847         PetscCall(PetscObjectTypeCompare((PetscObject)T, MATIS, &ismatis));
848         PetscCheck(ismatis, comm, PETSC_ERR_SUP, "Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") (transposed) is not of type MATIS", i, j);
849         PetscCall(MatISGetAllowRepeated(T, &sallow));
850         PetscCall(MatISGetLocalMat(T, &lT));
851         PetscCall(MatCreateTranspose(lT, &snest[ij]));
852       } else {
853         PetscCall(PetscObjectTypeCompare((PetscObject)nest[i][j], MATIS, &ismatis));
854         PetscCheck(ismatis, comm, PETSC_ERR_SUP, "Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") is not of type MATIS", i, j);
855         PetscCall(MatISGetAllowRepeated(nest[i][j], &sallow));
856         PetscCall(MatISGetLocalMat(nest[i][j], &snest[ij]));
857       }
858       if (allow_repeated == PETSC_BOOL3_UNKNOWN) allow_repeated = PetscBoolToBool3(sallow);
859       PetscCheck(sallow == PetscBool3ToBool(allow_repeated), comm, PETSC_ERR_SUP, "Cannot mix repeated and non repeated maps");
860 
861       /* Check compatibility of local sizes */
862       PetscCall(MatGetSize(snest[ij], &l1, &l2));
863       PetscCall(MatGetBlockSizes(snest[ij], &lb1, &lb2));
864       if (!l1 || !l2) continue;
865       PetscCheck(!lr[i] || l1 == lr[i], PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid local size %" PetscInt_FMT " != %" PetscInt_FMT, i, j, lr[i], l1);
866       PetscCheck(!lc[j] || l2 == lc[j], PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid local size %" PetscInt_FMT " != %" PetscInt_FMT, i, j, lc[j], l2);
867       lr[i] = l1;
868       lc[j] = l2;
869 
870       /* check compatibility for local matrix reusage */
871       if (rnest && !rnest[i][j] != !snest[ij]) lreuse = PETSC_FALSE;
872     }
873   }
874 
875   if (PetscDefined(USE_DEBUG)) {
876     /* Check compatibility of l2g maps for rows */
877     for (i = 0; i < nr; i++) {
878       rl2g = NULL;
879       for (j = 0; j < nc; j++) {
880         PetscInt n1, n2;
881 
882         if (!nest[i][j]) continue;
883         if (istrans[i * nc + j]) {
884           Mat T;
885 
886           PetscCall(MatTransposeGetMat(nest[i][j], &T));
887           PetscCall(MatISGetLocalToGlobalMapping(T, NULL, &cl2g));
888         } else {
889           PetscCall(MatISGetLocalToGlobalMapping(nest[i][j], &cl2g, NULL));
890         }
891         PetscCall(ISLocalToGlobalMappingGetSize(cl2g, &n1));
892         if (!n1) continue;
893         if (!rl2g) {
894           rl2g = cl2g;
895         } else {
896           const PetscInt *idxs1, *idxs2;
897           PetscBool       same;
898 
899           PetscCall(ISLocalToGlobalMappingGetSize(rl2g, &n2));
900           PetscCheck(n1 == n2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid row l2gmap size %" PetscInt_FMT " != %" PetscInt_FMT, i, j, n1, n2);
901           PetscCall(ISLocalToGlobalMappingGetIndices(cl2g, &idxs1));
902           PetscCall(ISLocalToGlobalMappingGetIndices(rl2g, &idxs2));
903           PetscCall(PetscArraycmp(idxs1, idxs2, n1, &same));
904           PetscCall(ISLocalToGlobalMappingRestoreIndices(cl2g, &idxs1));
905           PetscCall(ISLocalToGlobalMappingRestoreIndices(rl2g, &idxs2));
906           PetscCheck(same, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid row l2gmap", i, j);
907         }
908       }
909     }
910     /* Check compatibility of l2g maps for columns */
911     for (i = 0; i < nc; i++) {
912       rl2g = NULL;
913       for (j = 0; j < nr; j++) {
914         PetscInt n1, n2;
915 
916         if (!nest[j][i]) continue;
917         if (istrans[j * nc + i]) {
918           Mat T;
919 
920           PetscCall(MatTransposeGetMat(nest[j][i], &T));
921           PetscCall(MatISGetLocalToGlobalMapping(T, &cl2g, NULL));
922         } else {
923           PetscCall(MatISGetLocalToGlobalMapping(nest[j][i], NULL, &cl2g));
924         }
925         PetscCall(ISLocalToGlobalMappingGetSize(cl2g, &n1));
926         if (!n1) continue;
927         if (!rl2g) {
928           rl2g = cl2g;
929         } else {
930           const PetscInt *idxs1, *idxs2;
931           PetscBool       same;
932 
933           PetscCall(ISLocalToGlobalMappingGetSize(rl2g, &n2));
934           PetscCheck(n1 == n2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid column l2gmap size %" PetscInt_FMT " != %" PetscInt_FMT, j, i, n1, n2);
935           PetscCall(ISLocalToGlobalMappingGetIndices(cl2g, &idxs1));
936           PetscCall(ISLocalToGlobalMappingGetIndices(rl2g, &idxs2));
937           PetscCall(PetscArraycmp(idxs1, idxs2, n1, &same));
938           PetscCall(ISLocalToGlobalMappingRestoreIndices(cl2g, &idxs1));
939           PetscCall(ISLocalToGlobalMappingRestoreIndices(rl2g, &idxs2));
940           PetscCheck(same, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot convert from MATNEST to MATIS! Matrix block (%" PetscInt_FMT ",%" PetscInt_FMT ") has invalid column l2gmap", j, i);
941         }
942       }
943     }
944   }
945 
946   B = NULL;
947   if (reuse != MAT_REUSE_MATRIX) {
948     PetscInt stl;
949 
950     /* Create l2g map for the rows of the new matrix and index sets for the local MATNEST */
951     for (i = 0, stl = 0; i < nr; i++) stl += lr[i];
952     PetscCall(PetscMalloc1(stl, &l2gidxs));
953     for (i = 0, stl = 0; i < nr; i++) {
954       Mat             usedmat;
955       Mat_IS         *matis;
956       const PetscInt *idxs;
957 
958       /* local IS for local NEST */
959       PetscCall(ISCreateStride(PETSC_COMM_SELF, lr[i], stl, 1, &islrow[i]));
960 
961       /* l2gmap */
962       j       = 0;
963       usedmat = nest[i][j];
964       while (!usedmat && j < nc - 1) usedmat = nest[i][++j];
965       PetscCheck(usedmat, comm, PETSC_ERR_SUP, "Cannot find valid row mat");
966 
967       if (istrans[i * nc + j]) {
968         Mat T;
969         PetscCall(MatTransposeGetMat(usedmat, &T));
970         usedmat = T;
971       }
972       matis = (Mat_IS *)usedmat->data;
973       PetscCall(ISGetIndices(isrow[i], &idxs));
974       if (istrans[i * nc + j]) {
975         PetscCall(PetscSFBcastBegin(matis->csf, MPIU_INT, idxs, l2gidxs + stl, MPI_REPLACE));
976         PetscCall(PetscSFBcastEnd(matis->csf, MPIU_INT, idxs, l2gidxs + stl, MPI_REPLACE));
977       } else {
978         PetscCall(PetscSFBcastBegin(matis->sf, MPIU_INT, idxs, l2gidxs + stl, MPI_REPLACE));
979         PetscCall(PetscSFBcastEnd(matis->sf, MPIU_INT, idxs, l2gidxs + stl, MPI_REPLACE));
980       }
981       PetscCall(ISRestoreIndices(isrow[i], &idxs));
982       stl += lr[i];
983     }
984     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, stl, l2gidxs, PETSC_OWN_POINTER, &rl2g));
985 
986     /* Create l2g map for columns of the new matrix and index sets for the local MATNEST */
987     for (i = 0, stl = 0; i < nc; i++) stl += lc[i];
988     PetscCall(PetscMalloc1(stl, &l2gidxs));
989     for (i = 0, stl = 0; i < nc; i++) {
990       Mat             usedmat;
991       Mat_IS         *matis;
992       const PetscInt *idxs;
993 
994       /* local IS for local NEST */
995       PetscCall(ISCreateStride(PETSC_COMM_SELF, lc[i], stl, 1, &islcol[i]));
996 
997       /* l2gmap */
998       j       = 0;
999       usedmat = nest[j][i];
1000       while (!usedmat && j < nr - 1) usedmat = nest[++j][i];
1001       PetscCheck(usedmat, comm, PETSC_ERR_SUP, "Cannot find valid column mat");
1002       if (istrans[j * nc + i]) {
1003         Mat T;
1004         PetscCall(MatTransposeGetMat(usedmat, &T));
1005         usedmat = T;
1006       }
1007       matis = (Mat_IS *)usedmat->data;
1008       PetscCall(ISGetIndices(iscol[i], &idxs));
1009       if (istrans[j * nc + i]) {
1010         PetscCall(PetscSFBcastBegin(matis->sf, MPIU_INT, idxs, l2gidxs + stl, MPI_REPLACE));
1011         PetscCall(PetscSFBcastEnd(matis->sf, MPIU_INT, idxs, l2gidxs + stl, MPI_REPLACE));
1012       } else {
1013         PetscCall(PetscSFBcastBegin(matis->csf, MPIU_INT, idxs, l2gidxs + stl, MPI_REPLACE));
1014         PetscCall(PetscSFBcastEnd(matis->csf, MPIU_INT, idxs, l2gidxs + stl, MPI_REPLACE));
1015       }
1016       PetscCall(ISRestoreIndices(iscol[i], &idxs));
1017       stl += lc[i];
1018     }
1019     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, stl, l2gidxs, PETSC_OWN_POINTER, &cl2g));
1020 
1021     /* Create MATIS */
1022     PetscCall(MatCreate(comm, &B));
1023     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
1024     PetscCall(MatGetBlockSizes(A, &rbs, &cbs));
1025     PetscCall(MatSetBlockSizes(B, rbs, cbs));
1026     PetscCall(MatSetType(B, MATIS));
1027     PetscCall(MatISSetLocalMatType(B, MATNEST));
1028     PetscCall(MatISSetAllowRepeated(B, PetscBool3ToBool(allow_repeated)));
1029     { /* hack : avoid setup of scatters */
1030       Mat_IS *matis     = (Mat_IS *)B->data;
1031       matis->islocalref = PETSC_TRUE;
1032     }
1033     PetscCall(MatSetLocalToGlobalMapping(B, rl2g, cl2g));
1034     PetscCall(ISLocalToGlobalMappingDestroy(&rl2g));
1035     PetscCall(ISLocalToGlobalMappingDestroy(&cl2g));
1036     PetscCall(MatCreateNest(PETSC_COMM_SELF, nr, islrow, nc, islcol, snest, &lA));
1037     PetscCall(MatNestSetVecType(lA, VECNEST));
1038     for (i = 0; i < nr * nc; i++) {
1039       if (istrans[i]) PetscCall(MatDestroy(&snest[i]));
1040     }
1041     PetscCall(MatISSetLocalMat(B, lA));
1042     PetscCall(MatDestroy(&lA));
1043     { /* hack : setup of scatters done here */
1044       Mat_IS *matis = (Mat_IS *)B->data;
1045 
1046       matis->islocalref = PETSC_FALSE;
1047       PetscCall(MatISSetUpScatters_Private(B));
1048     }
1049     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1050     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1051     if (reuse == MAT_INPLACE_MATRIX) {
1052       PetscCall(MatHeaderReplace(A, &B));
1053     } else {
1054       *newmat = B;
1055     }
1056   } else {
1057     if (lreuse) {
1058       PetscCall(MatISGetLocalMat(*newmat, &lA));
1059       for (i = 0; i < nr; i++) {
1060         for (j = 0; j < nc; j++) {
1061           if (snest[i * nc + j]) {
1062             PetscCall(MatNestSetSubMat(lA, i, j, snest[i * nc + j]));
1063             if (istrans[i * nc + j]) PetscCall(MatDestroy(&snest[i * nc + j]));
1064           }
1065         }
1066       }
1067     } else {
1068       PetscInt stl;
1069       for (i = 0, stl = 0; i < nr; i++) {
1070         PetscCall(ISCreateStride(PETSC_COMM_SELF, lr[i], stl, 1, &islrow[i]));
1071         stl += lr[i];
1072       }
1073       for (i = 0, stl = 0; i < nc; i++) {
1074         PetscCall(ISCreateStride(PETSC_COMM_SELF, lc[i], stl, 1, &islcol[i]));
1075         stl += lc[i];
1076       }
1077       PetscCall(MatCreateNest(PETSC_COMM_SELF, nr, islrow, nc, islcol, snest, &lA));
1078       for (i = 0; i < nr * nc; i++) {
1079         if (istrans[i]) PetscCall(MatDestroy(&snest[i]));
1080       }
1081       PetscCall(MatISSetLocalMat(*newmat, lA));
1082       PetscCall(MatDestroy(&lA));
1083     }
1084     PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
1085     PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
1086   }
1087 
1088   /* Create local matrix in MATNEST format */
1089   convert = PETSC_FALSE;
1090   PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)A)->prefix, "-mat_is_convert_local_nest", &convert, NULL));
1091   if (convert) {
1092     Mat              M;
1093     MatISLocalFields lf;
1094 
1095     PetscCall(MatISGetLocalMat(*newmat, &lA));
1096     PetscCall(MatConvert(lA, MATAIJ, MAT_INITIAL_MATRIX, &M));
1097     PetscCall(MatISSetLocalMat(*newmat, M));
1098     PetscCall(MatDestroy(&M));
1099 
1100     /* attach local fields to the matrix */
1101     PetscCall(PetscNew(&lf));
1102     PetscCall(PetscMalloc2(nr, &lf->rf, nc, &lf->cf));
1103     for (i = 0; i < nr; i++) {
1104       PetscInt n, st;
1105 
1106       PetscCall(ISGetLocalSize(islrow[i], &n));
1107       PetscCall(ISStrideGetInfo(islrow[i], &st, NULL));
1108       PetscCall(ISCreateStride(comm, n, st, 1, &lf->rf[i]));
1109     }
1110     for (i = 0; i < nc; i++) {
1111       PetscInt n, st;
1112 
1113       PetscCall(ISGetLocalSize(islcol[i], &n));
1114       PetscCall(ISStrideGetInfo(islcol[i], &st, NULL));
1115       PetscCall(ISCreateStride(comm, n, st, 1, &lf->cf[i]));
1116     }
1117     lf->nr = nr;
1118     lf->nc = nc;
1119     PetscCall(PetscObjectContainerCompose((PetscObject)*newmat, "_convert_nest_lfields", lf, MatISContainerDestroyFields_Private));
1120   }
1121 
1122   /* Free workspace */
1123   for (i = 0; i < nr; i++) PetscCall(ISDestroy(&islrow[i]));
1124   for (i = 0; i < nc; i++) PetscCall(ISDestroy(&islcol[i]));
1125   PetscCall(PetscFree6(isrow, iscol, islrow, islcol, snest, istrans));
1126   PetscCall(PetscFree2(lr, lc));
1127   PetscFunctionReturn(PETSC_SUCCESS);
1128 }
1129 
1130 static PetscErrorCode MatDiagonalScale_IS(Mat A, Vec l, Vec r)
1131 {
1132   Mat_IS            *matis = (Mat_IS *)A->data;
1133   Vec                ll, rr;
1134   const PetscScalar *Y, *X;
1135   PetscScalar       *x, *y;
1136 
1137   PetscFunctionBegin;
1138   if (l) {
1139     ll = matis->y;
1140     PetscCall(VecGetArrayRead(l, &Y));
1141     PetscCall(VecGetArray(ll, &y));
1142     PetscCall(PetscSFBcastBegin(matis->sf, MPIU_SCALAR, Y, y, MPI_REPLACE));
1143   } else {
1144     ll = NULL;
1145   }
1146   if (r) {
1147     rr = matis->x;
1148     PetscCall(VecGetArrayRead(r, &X));
1149     PetscCall(VecGetArray(rr, &x));
1150     PetscCall(PetscSFBcastBegin(matis->csf, MPIU_SCALAR, X, x, MPI_REPLACE));
1151   } else {
1152     rr = NULL;
1153   }
1154   if (ll) {
1155     PetscCall(PetscSFBcastEnd(matis->sf, MPIU_SCALAR, Y, y, MPI_REPLACE));
1156     PetscCall(VecRestoreArrayRead(l, &Y));
1157     PetscCall(VecRestoreArray(ll, &y));
1158   }
1159   if (rr) {
1160     PetscCall(PetscSFBcastEnd(matis->csf, MPIU_SCALAR, X, x, MPI_REPLACE));
1161     PetscCall(VecRestoreArrayRead(r, &X));
1162     PetscCall(VecRestoreArray(rr, &x));
1163   }
1164   PetscCall(MatDiagonalScale(matis->A, ll, rr));
1165   PetscFunctionReturn(PETSC_SUCCESS);
1166 }
1167 
1168 static PetscErrorCode MatGetInfo_IS(Mat A, MatInfoType flag, MatInfo *ginfo)
1169 {
1170   Mat_IS        *matis = (Mat_IS *)A->data;
1171   MatInfo        info;
1172   PetscLogDouble isend[6], irecv[6];
1173   PetscInt       bs;
1174 
1175   PetscFunctionBegin;
1176   PetscCall(MatGetBlockSize(A, &bs));
1177   if (matis->A->ops->getinfo) {
1178     PetscCall(MatGetInfo(matis->A, MAT_LOCAL, &info));
1179     isend[0] = info.nz_used;
1180     isend[1] = info.nz_allocated;
1181     isend[2] = info.nz_unneeded;
1182     isend[3] = info.memory;
1183     isend[4] = info.mallocs;
1184   } else {
1185     isend[0] = 0.;
1186     isend[1] = 0.;
1187     isend[2] = 0.;
1188     isend[3] = 0.;
1189     isend[4] = 0.;
1190   }
1191   isend[5] = matis->A->num_ass;
1192   if (flag == MAT_LOCAL) {
1193     ginfo->nz_used      = isend[0];
1194     ginfo->nz_allocated = isend[1];
1195     ginfo->nz_unneeded  = isend[2];
1196     ginfo->memory       = isend[3];
1197     ginfo->mallocs      = isend[4];
1198     ginfo->assemblies   = isend[5];
1199   } else if (flag == MAT_GLOBAL_MAX) {
1200     PetscCallMPI(MPIU_Allreduce(isend, irecv, 6, MPIU_PETSCLOGDOUBLE, MPI_MAX, PetscObjectComm((PetscObject)A)));
1201 
1202     ginfo->nz_used      = irecv[0];
1203     ginfo->nz_allocated = irecv[1];
1204     ginfo->nz_unneeded  = irecv[2];
1205     ginfo->memory       = irecv[3];
1206     ginfo->mallocs      = irecv[4];
1207     ginfo->assemblies   = irecv[5];
1208   } else if (flag == MAT_GLOBAL_SUM) {
1209     PetscCallMPI(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_SUM, PetscObjectComm((PetscObject)A)));
1210 
1211     ginfo->nz_used      = irecv[0];
1212     ginfo->nz_allocated = irecv[1];
1213     ginfo->nz_unneeded  = irecv[2];
1214     ginfo->memory       = irecv[3];
1215     ginfo->mallocs      = irecv[4];
1216     ginfo->assemblies   = A->num_ass;
1217   }
1218   ginfo->block_size        = bs;
1219   ginfo->fill_ratio_given  = 0;
1220   ginfo->fill_ratio_needed = 0;
1221   ginfo->factor_mallocs    = 0;
1222   PetscFunctionReturn(PETSC_SUCCESS);
1223 }
1224 
1225 static PetscErrorCode MatTranspose_IS(Mat A, MatReuse reuse, Mat *B)
1226 {
1227   Mat C, lC, lA;
1228 
1229   PetscFunctionBegin;
1230   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
1231   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
1232     ISLocalToGlobalMapping rl2g, cl2g;
1233     PetscBool              allow_repeated;
1234 
1235     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
1236     PetscCall(MatSetSizes(C, A->cmap->n, A->rmap->n, A->cmap->N, A->rmap->N));
1237     PetscCall(MatSetBlockSizes(C, A->cmap->bs, A->rmap->bs));
1238     PetscCall(MatSetType(C, MATIS));
1239     PetscCall(MatISGetAllowRepeated(A, &allow_repeated));
1240     PetscCall(MatISSetAllowRepeated(C, allow_repeated));
1241     PetscCall(MatGetLocalToGlobalMapping(A, &rl2g, &cl2g));
1242     PetscCall(MatSetLocalToGlobalMapping(C, cl2g, rl2g));
1243   } else C = *B;
1244 
1245   /* perform local transposition */
1246   PetscCall(MatISGetLocalMat(A, &lA));
1247   PetscCall(MatTranspose(lA, MAT_INITIAL_MATRIX, &lC));
1248   PetscCall(MatSetLocalToGlobalMapping(lC, lA->cmap->mapping, lA->rmap->mapping));
1249   PetscCall(MatISSetLocalMat(C, lC));
1250   PetscCall(MatDestroy(&lC));
1251 
1252   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
1253     *B = C;
1254   } else {
1255     PetscCall(MatHeaderMerge(A, &C));
1256   }
1257   PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY));
1258   PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY));
1259   PetscFunctionReturn(PETSC_SUCCESS);
1260 }
1261 
1262 static PetscErrorCode MatDiagonalSet_IS(Mat A, Vec D, InsertMode insmode)
1263 {
1264   Mat_IS *is = (Mat_IS *)A->data;
1265 
1266   PetscFunctionBegin;
1267   PetscCheck(!is->allow_repeated || insmode == ADD_VALUES, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "INSERT_VALUES with repeated entries not supported");
1268   if (D) { /* MatShift_IS pass D = NULL */
1269     PetscCall(VecScatterBegin(is->rctx, D, is->y, INSERT_VALUES, SCATTER_FORWARD));
1270     PetscCall(VecScatterEnd(is->rctx, D, is->y, INSERT_VALUES, SCATTER_FORWARD));
1271   }
1272   PetscCall(VecPointwiseDivide(is->y, is->y, is->counter));
1273   PetscCall(MatDiagonalSet(is->A, is->y, insmode));
1274   PetscFunctionReturn(PETSC_SUCCESS);
1275 }
1276 
1277 static PetscErrorCode MatShift_IS(Mat A, PetscScalar a)
1278 {
1279   Mat_IS *is = (Mat_IS *)A->data;
1280 
1281   PetscFunctionBegin;
1282   PetscCall(VecSet(is->y, a));
1283   PetscCall(MatDiagonalSet_IS(A, NULL, ADD_VALUES));
1284   PetscFunctionReturn(PETSC_SUCCESS);
1285 }
1286 
1287 static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode addv)
1288 {
1289   PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION], cols_l[MATIS_MAX_ENTRIES_INSERTION];
1290 
1291   PetscFunctionBegin;
1292   PetscCheck(m <= MATIS_MAX_ENTRIES_INSERTION && n <= MATIS_MAX_ENTRIES_INSERTION, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of row/column indices must be <= %d: they are %" PetscInt_FMT " %" PetscInt_FMT, MATIS_MAX_ENTRIES_INSERTION, m, n);
1293   PetscCall(ISLocalToGlobalMappingApply(A->rmap->mapping, m, rows, rows_l));
1294   PetscCall(ISLocalToGlobalMappingApply(A->cmap->mapping, n, cols, cols_l));
1295   PetscCall(MatSetValuesLocal_IS(A, m, rows_l, n, cols_l, values, addv));
1296   PetscFunctionReturn(PETSC_SUCCESS);
1297 }
1298 
1299 static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode addv)
1300 {
1301   PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION], cols_l[MATIS_MAX_ENTRIES_INSERTION];
1302 
1303   PetscFunctionBegin;
1304   PetscCheck(m <= MATIS_MAX_ENTRIES_INSERTION && n <= MATIS_MAX_ENTRIES_INSERTION, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of row/column block indices must be <= %d: they are %" PetscInt_FMT " %" PetscInt_FMT, MATIS_MAX_ENTRIES_INSERTION, m, n);
1305   PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, m, rows, rows_l));
1306   PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, n, cols, cols_l));
1307   PetscCall(MatSetValuesBlockedLocal_IS(A, m, rows_l, n, cols_l, values, addv));
1308   PetscFunctionReturn(PETSC_SUCCESS);
1309 }
1310 
1311 static PetscErrorCode MatCreateSubMatrix_IS(Mat mat, IS irow, IS icol, MatReuse scall, Mat *newmat)
1312 {
1313   Mat             locmat, newlocmat;
1314   Mat_IS         *newmatis;
1315   const PetscInt *idxs;
1316   PetscInt        i, m, n;
1317 
1318   PetscFunctionBegin;
1319   if (scall == MAT_REUSE_MATRIX) {
1320     PetscBool ismatis;
1321 
1322     PetscCall(PetscObjectTypeCompare((PetscObject)*newmat, MATIS, &ismatis));
1323     PetscCheck(ismatis, PetscObjectComm((PetscObject)*newmat), PETSC_ERR_ARG_WRONG, "Cannot reuse matrix! Not of MATIS type");
1324     newmatis = (Mat_IS *)(*newmat)->data;
1325     PetscCheck(newmatis->getsub_ris, PetscObjectComm((PetscObject)*newmat), PETSC_ERR_ARG_WRONG, "Cannot reuse matrix! Misses local row IS");
1326     PetscCheck(newmatis->getsub_cis, PetscObjectComm((PetscObject)*newmat), PETSC_ERR_ARG_WRONG, "Cannot reuse matrix! Misses local col IS");
1327   }
1328   /* irow and icol may not have duplicate entries */
1329   if (PetscDefined(USE_DEBUG)) {
1330     Vec                rtest, ltest;
1331     const PetscScalar *array;
1332 
1333     PetscCall(MatCreateVecs(mat, &ltest, &rtest));
1334     PetscCall(ISGetLocalSize(irow, &n));
1335     PetscCall(ISGetIndices(irow, &idxs));
1336     for (i = 0; i < n; i++) PetscCall(VecSetValue(rtest, idxs[i], 1.0, ADD_VALUES));
1337     PetscCall(VecAssemblyBegin(rtest));
1338     PetscCall(VecAssemblyEnd(rtest));
1339     PetscCall(VecGetLocalSize(rtest, &n));
1340     PetscCall(VecGetOwnershipRange(rtest, &m, NULL));
1341     PetscCall(VecGetArrayRead(rtest, &array));
1342     for (i = 0; i < n; i++) PetscCheck(array[i] == 0. || array[i] == 1., PETSC_COMM_SELF, PETSC_ERR_SUP, "Index %" PetscInt_FMT " counted %" PetscInt_FMT " times! Irow may not have duplicate entries", i + m, (PetscInt)PetscRealPart(array[i]));
1343     PetscCall(VecRestoreArrayRead(rtest, &array));
1344     PetscCall(ISRestoreIndices(irow, &idxs));
1345     PetscCall(ISGetLocalSize(icol, &n));
1346     PetscCall(ISGetIndices(icol, &idxs));
1347     for (i = 0; i < n; i++) PetscCall(VecSetValue(ltest, idxs[i], 1.0, ADD_VALUES));
1348     PetscCall(VecAssemblyBegin(ltest));
1349     PetscCall(VecAssemblyEnd(ltest));
1350     PetscCall(VecGetLocalSize(ltest, &n));
1351     PetscCall(VecGetOwnershipRange(ltest, &m, NULL));
1352     PetscCall(VecGetArrayRead(ltest, &array));
1353     for (i = 0; i < n; i++) PetscCheck(array[i] == 0. || array[i] == 1., PETSC_COMM_SELF, PETSC_ERR_SUP, "Index %" PetscInt_FMT " counted %" PetscInt_FMT " times! Icol may not have duplicate entries", i + m, (PetscInt)PetscRealPart(array[i]));
1354     PetscCall(VecRestoreArrayRead(ltest, &array));
1355     PetscCall(ISRestoreIndices(icol, &idxs));
1356     PetscCall(VecDestroy(&rtest));
1357     PetscCall(VecDestroy(&ltest));
1358   }
1359   if (scall == MAT_INITIAL_MATRIX) {
1360     Mat_IS                *matis = (Mat_IS *)mat->data;
1361     ISLocalToGlobalMapping rl2g;
1362     IS                     is;
1363     PetscInt              *lidxs, *lgidxs, *newgidxs;
1364     PetscInt               ll, newloc, irbs, icbs, arbs, acbs, rbs, cbs;
1365     PetscBool              cong;
1366     MPI_Comm               comm;
1367 
1368     PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
1369     PetscCall(MatGetBlockSizes(mat, &arbs, &acbs));
1370     PetscCall(ISGetBlockSize(irow, &irbs));
1371     PetscCall(ISGetBlockSize(icol, &icbs));
1372     rbs = arbs == irbs ? irbs : 1;
1373     cbs = acbs == icbs ? icbs : 1;
1374     PetscCall(ISGetLocalSize(irow, &m));
1375     PetscCall(ISGetLocalSize(icol, &n));
1376     PetscCall(MatCreate(comm, newmat));
1377     PetscCall(MatSetType(*newmat, MATIS));
1378     PetscCall(MatISSetAllowRepeated(*newmat, matis->allow_repeated));
1379     PetscCall(MatSetSizes(*newmat, m, n, PETSC_DECIDE, PETSC_DECIDE));
1380     PetscCall(MatSetBlockSizes(*newmat, rbs, cbs));
1381     /* communicate irow to their owners in the layout */
1382     PetscCall(ISGetIndices(irow, &idxs));
1383     PetscCall(PetscLayoutMapLocal(mat->rmap, m, idxs, &ll, &lidxs, &lgidxs));
1384     PetscCall(ISRestoreIndices(irow, &idxs));
1385     PetscCall(PetscArrayzero(matis->sf_rootdata, matis->sf->nroots));
1386     for (i = 0; i < ll; i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i] + 1;
1387     PetscCall(PetscFree(lidxs));
1388     PetscCall(PetscFree(lgidxs));
1389     PetscCall(PetscSFBcastBegin(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
1390     PetscCall(PetscSFBcastEnd(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
1391     for (i = 0, newloc = 0; i < matis->sf->nleaves; i++)
1392       if (matis->sf_leafdata[i]) newloc++;
1393     PetscCall(PetscMalloc1(newloc, &newgidxs));
1394     PetscCall(PetscMalloc1(newloc, &lidxs));
1395     for (i = 0, newloc = 0; i < matis->sf->nleaves; i++)
1396       if (matis->sf_leafdata[i]) {
1397         lidxs[newloc]      = i;
1398         newgidxs[newloc++] = matis->sf_leafdata[i] - 1;
1399       }
1400     PetscCall(ISCreateGeneral(comm, newloc, newgidxs, PETSC_OWN_POINTER, &is));
1401     PetscCall(ISLocalToGlobalMappingCreateIS(is, &rl2g));
1402     PetscCall(ISLocalToGlobalMappingSetBlockSize(rl2g, rbs));
1403     PetscCall(ISDestroy(&is));
1404     /* local is to extract local submatrix */
1405     newmatis = (Mat_IS *)(*newmat)->data;
1406     PetscCall(ISCreateGeneral(comm, newloc, lidxs, PETSC_OWN_POINTER, &newmatis->getsub_ris));
1407     PetscCall(MatHasCongruentLayouts(mat, &cong));
1408     if (cong && irow == icol && matis->csf == matis->sf) {
1409       PetscCall(MatSetLocalToGlobalMapping(*newmat, rl2g, rl2g));
1410       PetscCall(PetscObjectReference((PetscObject)newmatis->getsub_ris));
1411       newmatis->getsub_cis = newmatis->getsub_ris;
1412     } else {
1413       ISLocalToGlobalMapping cl2g;
1414 
1415       /* communicate icol to their owners in the layout */
1416       PetscCall(ISGetIndices(icol, &idxs));
1417       PetscCall(PetscLayoutMapLocal(mat->cmap, n, idxs, &ll, &lidxs, &lgidxs));
1418       PetscCall(ISRestoreIndices(icol, &idxs));
1419       PetscCall(PetscArrayzero(matis->csf_rootdata, matis->csf->nroots));
1420       for (i = 0; i < ll; i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i] + 1;
1421       PetscCall(PetscFree(lidxs));
1422       PetscCall(PetscFree(lgidxs));
1423       PetscCall(PetscSFBcastBegin(matis->csf, MPIU_INT, matis->csf_rootdata, matis->csf_leafdata, MPI_REPLACE));
1424       PetscCall(PetscSFBcastEnd(matis->csf, MPIU_INT, matis->csf_rootdata, matis->csf_leafdata, MPI_REPLACE));
1425       for (i = 0, newloc = 0; i < matis->csf->nleaves; i++)
1426         if (matis->csf_leafdata[i]) newloc++;
1427       PetscCall(PetscMalloc1(newloc, &newgidxs));
1428       PetscCall(PetscMalloc1(newloc, &lidxs));
1429       for (i = 0, newloc = 0; i < matis->csf->nleaves; i++)
1430         if (matis->csf_leafdata[i]) {
1431           lidxs[newloc]      = i;
1432           newgidxs[newloc++] = matis->csf_leafdata[i] - 1;
1433         }
1434       PetscCall(ISCreateGeneral(comm, newloc, newgidxs, PETSC_OWN_POINTER, &is));
1435       PetscCall(ISLocalToGlobalMappingCreateIS(is, &cl2g));
1436       PetscCall(ISLocalToGlobalMappingSetBlockSize(cl2g, cbs));
1437       PetscCall(ISDestroy(&is));
1438       /* local is to extract local submatrix */
1439       PetscCall(ISCreateGeneral(comm, newloc, lidxs, PETSC_OWN_POINTER, &newmatis->getsub_cis));
1440       PetscCall(MatSetLocalToGlobalMapping(*newmat, rl2g, cl2g));
1441       PetscCall(ISLocalToGlobalMappingDestroy(&cl2g));
1442     }
1443     PetscCall(ISLocalToGlobalMappingDestroy(&rl2g));
1444   } else {
1445     PetscCall(MatISGetLocalMat(*newmat, &newlocmat));
1446   }
1447   PetscCall(MatISGetLocalMat(mat, &locmat));
1448   newmatis = (Mat_IS *)(*newmat)->data;
1449   PetscCall(MatCreateSubMatrix(locmat, newmatis->getsub_ris, newmatis->getsub_cis, scall, &newlocmat));
1450   if (scall == MAT_INITIAL_MATRIX) {
1451     PetscCall(MatISSetLocalMat(*newmat, newlocmat));
1452     PetscCall(MatDestroy(&newlocmat));
1453   }
1454   PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
1455   PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
1456   PetscFunctionReturn(PETSC_SUCCESS);
1457 }
1458 
1459 static PetscErrorCode MatCopy_IS(Mat A, Mat B, MatStructure str)
1460 {
1461   Mat_IS   *a = (Mat_IS *)A->data, *b;
1462   PetscBool ismatis;
1463 
1464   PetscFunctionBegin;
1465   PetscCall(PetscObjectTypeCompare((PetscObject)B, MATIS, &ismatis));
1466   PetscCheck(ismatis, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "Need to be implemented");
1467   b = (Mat_IS *)B->data;
1468   PetscCall(MatCopy(a->A, b->A, str));
1469   PetscCall(PetscObjectStateIncrease((PetscObject)B));
1470   PetscFunctionReturn(PETSC_SUCCESS);
1471 }
1472 
1473 static PetscErrorCode MatMissingDiagonal_IS(Mat A, PetscBool *missing, PetscInt *d)
1474 {
1475   Vec                v;
1476   const PetscScalar *array;
1477   PetscInt           i, n;
1478 
1479   PetscFunctionBegin;
1480   *missing = PETSC_FALSE;
1481   PetscCall(MatCreateVecs(A, NULL, &v));
1482   PetscCall(MatGetDiagonal(A, v));
1483   PetscCall(VecGetLocalSize(v, &n));
1484   PetscCall(VecGetArrayRead(v, &array));
1485   for (i = 0; i < n; i++)
1486     if (array[i] == 0.) break;
1487   PetscCall(VecRestoreArrayRead(v, &array));
1488   PetscCall(VecDestroy(&v));
1489   if (i != n) *missing = PETSC_TRUE;
1490   if (d) {
1491     *d = -1;
1492     if (*missing) {
1493       PetscInt rstart;
1494       PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
1495       *d = i + rstart;
1496     }
1497   }
1498   PetscFunctionReturn(PETSC_SUCCESS);
1499 }
1500 
1501 static PetscErrorCode MatISSetUpSF_IS(Mat B)
1502 {
1503   Mat_IS         *matis = (Mat_IS *)B->data;
1504   const PetscInt *gidxs;
1505   PetscInt        nleaves;
1506 
1507   PetscFunctionBegin;
1508   if (matis->sf) PetscFunctionReturn(PETSC_SUCCESS);
1509   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)B), &matis->sf));
1510   PetscCall(ISLocalToGlobalMappingGetIndices(matis->rmapping, &gidxs));
1511   PetscCall(ISLocalToGlobalMappingGetSize(matis->rmapping, &nleaves));
1512   PetscCall(PetscSFSetGraphLayout(matis->sf, B->rmap, nleaves, NULL, PETSC_OWN_POINTER, gidxs));
1513   PetscCall(ISLocalToGlobalMappingRestoreIndices(matis->rmapping, &gidxs));
1514   PetscCall(PetscMalloc2(matis->sf->nroots, &matis->sf_rootdata, matis->sf->nleaves, &matis->sf_leafdata));
1515   if (matis->rmapping != matis->cmapping) { /* setup SF for columns */
1516     PetscCall(ISLocalToGlobalMappingGetSize(matis->cmapping, &nleaves));
1517     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)B), &matis->csf));
1518     PetscCall(ISLocalToGlobalMappingGetIndices(matis->cmapping, &gidxs));
1519     PetscCall(PetscSFSetGraphLayout(matis->csf, B->cmap, nleaves, NULL, PETSC_OWN_POINTER, gidxs));
1520     PetscCall(ISLocalToGlobalMappingRestoreIndices(matis->cmapping, &gidxs));
1521     PetscCall(PetscMalloc2(matis->csf->nroots, &matis->csf_rootdata, matis->csf->nleaves, &matis->csf_leafdata));
1522   } else {
1523     matis->csf          = matis->sf;
1524     matis->csf_leafdata = matis->sf_leafdata;
1525     matis->csf_rootdata = matis->sf_rootdata;
1526   }
1527   PetscFunctionReturn(PETSC_SUCCESS);
1528 }
1529 
1530 /*@
1531   MatISGetAllowRepeated - Get the flag to allow repeated entries in the local to global map
1532 
1533   Not Collective
1534 
1535   Input Parameter:
1536 . A - the matrix
1537 
1538   Output Parameter:
1539 . flg - the boolean flag
1540 
1541   Level: intermediate
1542 
1543 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateIS()`, `MatSetLocalToGlobalMapping()`, `MatISSetAllowRepeated()`
1544 @*/
1545 PetscErrorCode MatISGetAllowRepeated(Mat A, PetscBool *flg)
1546 {
1547   PetscBool ismatis;
1548 
1549   PetscFunctionBegin;
1550   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1551   PetscAssertPointer(flg, 2);
1552   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATIS, &ismatis));
1553   PetscCheck(ismatis, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for matrix type %s", ((PetscObject)A)->type_name);
1554   *flg = ((Mat_IS *)A->data)->allow_repeated;
1555   PetscFunctionReturn(PETSC_SUCCESS);
1556 }
1557 
1558 /*@
1559   MatISSetAllowRepeated - Set the flag to allow repeated entries in the local to global map
1560 
1561   Logically Collective
1562 
1563   Input Parameters:
1564 + A   - the matrix
1565 - flg - the boolean flag
1566 
1567   Level: intermediate
1568 
1569   Notes:
1570   The default value is `PETSC_FALSE`.
1571   When called AFTER calling `MatSetLocalToGlobalMapping()` it will recreate the local matrices
1572   if `flg` is different from the previously set value.
1573   Specifically, when `flg` is true it will just recreate the local matrices, while if
1574   `flg` is false will assemble the local matrices summing up repeated entries.
1575 
1576 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateIS()`, `MatSetLocalToGlobalMapping()`, `MatISGetAllowRepeated()`
1577 @*/
1578 PetscErrorCode MatISSetAllowRepeated(Mat A, PetscBool flg)
1579 {
1580   PetscFunctionBegin;
1581   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1582   PetscValidType(A, 1);
1583   PetscValidLogicalCollectiveBool(A, flg, 2);
1584   PetscTryMethod(A, "MatISSetAllowRepeated_C", (Mat, PetscBool), (A, flg));
1585   PetscFunctionReturn(PETSC_SUCCESS);
1586 }
1587 
1588 static PetscErrorCode MatISSetAllowRepeated_IS(Mat A, PetscBool flg)
1589 {
1590   Mat_IS                *matis = (Mat_IS *)A->data;
1591   Mat                    lA    = NULL;
1592   ISLocalToGlobalMapping lrmap, lcmap;
1593 
1594   PetscFunctionBegin;
1595   if (flg == matis->allow_repeated) PetscFunctionReturn(PETSC_SUCCESS);
1596   if (!matis->A) { /* matrix has not been preallocated yet */
1597     matis->allow_repeated = flg;
1598     PetscFunctionReturn(PETSC_SUCCESS);
1599   }
1600   PetscCheck(!matis->islocalref, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented for local references");
1601   if (matis->allow_repeated) { /* we will assemble the old local matrix if needed */
1602     lA = matis->A;
1603     PetscCall(PetscObjectReference((PetscObject)lA));
1604   }
1605   /* In case flg is True, we only recreate the local matrix */
1606   matis->allow_repeated = flg;
1607   PetscCall(MatSetLocalToGlobalMapping(A, A->rmap->mapping, A->cmap->mapping));
1608   if (lA) { /* assemble previous local matrix if needed */
1609     Mat nA = matis->A;
1610 
1611     PetscCall(MatGetLocalToGlobalMapping(nA, &lrmap, &lcmap));
1612     if (!lrmap && !lcmap) {
1613       PetscCall(MatISSetLocalMat(A, lA));
1614     } else {
1615       Mat            P = NULL, R = NULL;
1616       MatProductType ptype;
1617 
1618       if (lrmap == lcmap) {
1619         ptype = MATPRODUCT_PtAP;
1620         PetscCall(MatCreateFromISLocalToGlobalMapping(lcmap, nA, PETSC_TRUE, PETSC_FALSE, NULL, &P));
1621         PetscCall(MatProductCreate(lA, P, NULL, &nA));
1622       } else {
1623         if (lcmap) PetscCall(MatCreateFromISLocalToGlobalMapping(lcmap, nA, PETSC_TRUE, PETSC_FALSE, NULL, &P));
1624         if (lrmap) PetscCall(MatCreateFromISLocalToGlobalMapping(lrmap, nA, PETSC_FALSE, PETSC_TRUE, NULL, &R));
1625         if (R && P) {
1626           ptype = MATPRODUCT_ABC;
1627           PetscCall(MatProductCreate(R, lA, P, &nA));
1628         } else if (R) {
1629           ptype = MATPRODUCT_AB;
1630           PetscCall(MatProductCreate(R, lA, NULL, &nA));
1631         } else {
1632           ptype = MATPRODUCT_AB;
1633           PetscCall(MatProductCreate(lA, P, NULL, &nA));
1634         }
1635       }
1636       PetscCall(MatProductSetType(nA, ptype));
1637       PetscCall(MatProductSetFromOptions(nA));
1638       PetscCall(MatProductSymbolic(nA));
1639       PetscCall(MatProductNumeric(nA));
1640       PetscCall(MatProductClear(nA));
1641       PetscCall(MatConvert(nA, matis->lmattype, MAT_INPLACE_MATRIX, &nA));
1642       PetscCall(MatISSetLocalMat(A, nA));
1643       PetscCall(MatDestroy(&nA));
1644       PetscCall(MatDestroy(&P));
1645       PetscCall(MatDestroy(&R));
1646     }
1647   }
1648   PetscCall(MatDestroy(&lA));
1649   PetscFunctionReturn(PETSC_SUCCESS);
1650 }
1651 
1652 /*@
1653   MatISStoreL2L - Store local-to-local operators during the Galerkin process of computing `MatPtAP()`
1654 
1655   Logically Collective
1656 
1657   Input Parameters:
1658 + A     - the matrix
1659 - store - the boolean flag
1660 
1661   Level: advanced
1662 
1663 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateIS()`, `MatISSetPreallocation()`, `MatPtAP()`
1664 @*/
1665 PetscErrorCode MatISStoreL2L(Mat A, PetscBool store)
1666 {
1667   PetscFunctionBegin;
1668   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1669   PetscValidType(A, 1);
1670   PetscValidLogicalCollectiveBool(A, store, 2);
1671   PetscTryMethod(A, "MatISStoreL2L_C", (Mat, PetscBool), (A, store));
1672   PetscFunctionReturn(PETSC_SUCCESS);
1673 }
1674 
1675 static PetscErrorCode MatISStoreL2L_IS(Mat A, PetscBool store)
1676 {
1677   Mat_IS *matis = (Mat_IS *)A->data;
1678 
1679   PetscFunctionBegin;
1680   matis->storel2l = store;
1681   if (!store) PetscCall(PetscObjectCompose((PetscObject)A, "_MatIS_PtAP_l2l", NULL));
1682   PetscFunctionReturn(PETSC_SUCCESS);
1683 }
1684 
1685 /*@
1686   MatISFixLocalEmpty - Compress out zero local rows from the local matrices
1687 
1688   Logically Collective
1689 
1690   Input Parameters:
1691 + A   - the matrix
1692 - fix - the boolean flag
1693 
1694   Level: advanced
1695 
1696   Note:
1697   When `fix` is `PETSC_TRUE`, new local matrices and l2g maps are generated during the final assembly process.
1698 
1699 .seealso: [](ch_matrices), `Mat`, `MATIS`, `MatCreate()`, `MatCreateIS()`, `MatISSetPreallocation()`, `MatAssemblyEnd()`, `MAT_FINAL_ASSEMBLY`
1700 @*/
1701 PetscErrorCode MatISFixLocalEmpty(Mat A, PetscBool fix)
1702 {
1703   PetscFunctionBegin;
1704   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1705   PetscValidType(A, 1);
1706   PetscValidLogicalCollectiveBool(A, fix, 2);
1707   PetscTryMethod(A, "MatISFixLocalEmpty_C", (Mat, PetscBool), (A, fix));
1708   PetscFunctionReturn(PETSC_SUCCESS);
1709 }
1710 
1711 static PetscErrorCode MatISFixLocalEmpty_IS(Mat A, PetscBool fix)
1712 {
1713   Mat_IS *matis = (Mat_IS *)A->data;
1714 
1715   PetscFunctionBegin;
1716   matis->locempty = fix;
1717   PetscFunctionReturn(PETSC_SUCCESS);
1718 }
1719 
1720 /*@
1721   MatISSetPreallocation - Preallocates memory for a `MATIS` parallel matrix.
1722 
1723   Collective
1724 
1725   Input Parameters:
1726 + B     - the matrix
1727 . d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1728            (same value is used for all local rows)
1729 . d_nnz - array containing the number of nonzeros in the various rows of the
1730            DIAGONAL portion of the local submatrix (possibly different for each row)
1731            or `NULL`, if `d_nz` is used to specify the nonzero structure.
1732            The size of this array is equal to the number of local rows, i.e `m`.
1733            For matrices that will be factored, you must leave room for (and set)
1734            the diagonal entry even if it is zero.
1735 . o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1736            submatrix (same value is used for all local rows).
1737 - o_nnz - array containing the number of nonzeros in the various rows of the
1738            OFF-DIAGONAL portion of the local submatrix (possibly different for
1739            each row) or `NULL`, if `o_nz` is used to specify the nonzero
1740            structure. The size of this array is equal to the number
1741            of local rows, i.e `m`.
1742 
1743    If the *_nnz parameter is given then the *_nz parameter is ignored
1744 
1745   Level: intermediate
1746 
1747   Note:
1748   This function has the same interface as the `MATMPIAIJ` preallocation routine in order to simplify the transition
1749   from the asssembled format to the unassembled one. It overestimates the preallocation of `MATIS` local
1750   matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
1751 
1752 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateIS()`, `MatMPIAIJSetPreallocation()`, `MatISGetLocalMat()`, `MATIS`
1753 @*/
1754 PetscErrorCode MatISSetPreallocation(Mat B, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[])
1755 {
1756   PetscFunctionBegin;
1757   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
1758   PetscValidType(B, 1);
1759   PetscTryMethod(B, "MatISSetPreallocation_C", (Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]), (B, d_nz, d_nnz, o_nz, o_nnz));
1760   PetscFunctionReturn(PETSC_SUCCESS);
1761 }
1762 
1763 static PetscErrorCode MatISSetPreallocation_IS(Mat B, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[])
1764 {
1765   Mat_IS  *matis = (Mat_IS *)B->data;
1766   PetscInt bs, i, nlocalcols;
1767 
1768   PetscFunctionBegin;
1769   PetscCall(MatSetUp(B));
1770   if (!d_nnz)
1771     for (i = 0; i < matis->sf->nroots; i++) matis->sf_rootdata[i] = d_nz;
1772   else
1773     for (i = 0; i < matis->sf->nroots; i++) matis->sf_rootdata[i] = d_nnz[i];
1774 
1775   if (!o_nnz)
1776     for (i = 0; i < matis->sf->nroots; i++) matis->sf_rootdata[i] += o_nz;
1777   else
1778     for (i = 0; i < matis->sf->nroots; i++) matis->sf_rootdata[i] += o_nnz[i];
1779 
1780   PetscCall(PetscSFBcastBegin(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
1781   PetscCall(MatGetSize(matis->A, NULL, &nlocalcols));
1782   PetscCall(MatGetBlockSize(matis->A, &bs));
1783   PetscCall(PetscSFBcastEnd(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
1784 
1785   for (i = 0; i < matis->sf->nleaves; i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i], nlocalcols);
1786   PetscCall(MatSeqAIJSetPreallocation(matis->A, 0, matis->sf_leafdata));
1787 #if defined(PETSC_HAVE_HYPRE)
1788   PetscCall(MatHYPRESetPreallocation(matis->A, 0, matis->sf_leafdata, 0, NULL));
1789 #endif
1790 
1791   for (i = 0; i < matis->sf->nleaves / bs; i++) {
1792     PetscInt b;
1793 
1794     matis->sf_leafdata[i] = matis->sf_leafdata[i * bs] / bs;
1795     for (b = 1; b < bs; b++) matis->sf_leafdata[i] = PetscMax(matis->sf_leafdata[i], matis->sf_leafdata[i * bs + b] / bs);
1796   }
1797   PetscCall(MatSeqBAIJSetPreallocation(matis->A, bs, 0, matis->sf_leafdata));
1798 
1799   nlocalcols /= bs;
1800   for (i = 0; i < matis->sf->nleaves / bs; i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i], nlocalcols - i);
1801   PetscCall(MatSeqSBAIJSetPreallocation(matis->A, bs, 0, matis->sf_leafdata));
1802 
1803   /* for other matrix types */
1804   PetscCall(MatSetUp(matis->A));
1805   PetscFunctionReturn(PETSC_SUCCESS);
1806 }
1807 
1808 PETSC_SINGLE_LIBRARY_VISIBILITY_INTERNAL PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
1809 {
1810   Mat_IS         *matis = (Mat_IS *)A->data;
1811   PetscInt       *my_dnz, *my_onz, *dnz, *onz, *mat_ranges, *row_ownership;
1812   const PetscInt *global_indices_r, *global_indices_c;
1813   PetscInt        i, j, bs, rows, cols;
1814   PetscInt        lrows, lcols;
1815   PetscInt        local_rows, local_cols;
1816   PetscMPIInt     size;
1817   PetscBool       isdense, issbaij;
1818 
1819   PetscFunctionBegin;
1820   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1821   PetscCall(MatGetSize(A, &rows, &cols));
1822   PetscCall(MatGetBlockSize(A, &bs));
1823   PetscCall(MatGetSize(matis->A, &local_rows, &local_cols));
1824   PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A, MATSEQDENSE, &isdense));
1825   PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A, MATSEQSBAIJ, &issbaij));
1826   PetscCall(ISLocalToGlobalMappingGetIndices(matis->rmapping, &global_indices_r));
1827   if (matis->rmapping != matis->cmapping) {
1828     PetscCall(ISLocalToGlobalMappingGetIndices(matis->cmapping, &global_indices_c));
1829   } else global_indices_c = global_indices_r;
1830 
1831   if (issbaij) PetscCall(MatGetRowUpperTriangular(matis->A));
1832   /*
1833      An SF reduce is needed to sum up properly on shared rows.
1834      Note that generally preallocation is not exact, since it overestimates nonzeros
1835   */
1836   PetscCall(MatGetLocalSize(A, &lrows, &lcols));
1837   MatPreallocateBegin(PetscObjectComm((PetscObject)A), lrows, lcols, dnz, onz);
1838   /* All processes need to compute entire row ownership */
1839   PetscCall(PetscMalloc1(rows, &row_ownership));
1840   PetscCall(MatGetOwnershipRanges(A, (const PetscInt **)&mat_ranges));
1841   for (i = 0; i < size; i++) {
1842     for (j = mat_ranges[i]; j < mat_ranges[i + 1]; j++) row_ownership[j] = i;
1843   }
1844   PetscCall(MatGetOwnershipRangesColumn(A, (const PetscInt **)&mat_ranges));
1845 
1846   /*
1847      my_dnz and my_onz contains exact contribution to preallocation from each local mat
1848      then, they will be summed up properly. This way, preallocation is always sufficient
1849   */
1850   PetscCall(PetscCalloc2(local_rows, &my_dnz, local_rows, &my_onz));
1851   /* preallocation as a MATAIJ */
1852   if (isdense) { /* special case for dense local matrices */
1853     for (i = 0; i < local_rows; i++) {
1854       PetscInt owner = row_ownership[global_indices_r[i]];
1855       for (j = 0; j < local_cols; j++) {
1856         PetscInt index_col = global_indices_c[j];
1857         if (index_col > mat_ranges[owner] - 1 && index_col < mat_ranges[owner + 1]) { /* diag block */
1858           my_dnz[i] += 1;
1859         } else { /* offdiag block */
1860           my_onz[i] += 1;
1861         }
1862       }
1863     }
1864   } else if (matis->A->ops->getrowij) {
1865     const PetscInt *ii, *jj, *jptr;
1866     PetscBool       done;
1867     PetscCall(MatGetRowIJ(matis->A, 0, PETSC_FALSE, PETSC_FALSE, &local_rows, &ii, &jj, &done));
1868     PetscCheck(done, PetscObjectComm((PetscObject)matis->A), PETSC_ERR_PLIB, "Error in MatGetRowIJ");
1869     jptr = jj;
1870     for (i = 0; i < local_rows; i++) {
1871       PetscInt index_row = global_indices_r[i];
1872       for (j = 0; j < ii[i + 1] - ii[i]; j++, jptr++) {
1873         PetscInt owner     = row_ownership[index_row];
1874         PetscInt index_col = global_indices_c[*jptr];
1875         if (index_col > mat_ranges[owner] - 1 && index_col < mat_ranges[owner + 1]) { /* diag block */
1876           my_dnz[i] += 1;
1877         } else { /* offdiag block */
1878           my_onz[i] += 1;
1879         }
1880         /* same as before, interchanging rows and cols */
1881         if (issbaij && index_col != index_row) {
1882           owner = row_ownership[index_col];
1883           if (index_row > mat_ranges[owner] - 1 && index_row < mat_ranges[owner + 1]) {
1884             my_dnz[*jptr] += 1;
1885           } else {
1886             my_onz[*jptr] += 1;
1887           }
1888         }
1889       }
1890     }
1891     PetscCall(MatRestoreRowIJ(matis->A, 0, PETSC_FALSE, PETSC_FALSE, &local_rows, &ii, &jj, &done));
1892     PetscCheck(done, PetscObjectComm((PetscObject)matis->A), PETSC_ERR_PLIB, "Error in MatRestoreRowIJ");
1893   } else { /* loop over rows and use MatGetRow */
1894     for (i = 0; i < local_rows; i++) {
1895       const PetscInt *cols;
1896       PetscInt        ncols, index_row = global_indices_r[i];
1897       PetscCall(MatGetRow(matis->A, i, &ncols, &cols, NULL));
1898       for (j = 0; j < ncols; j++) {
1899         PetscInt owner     = row_ownership[index_row];
1900         PetscInt index_col = global_indices_c[cols[j]];
1901         if (index_col > mat_ranges[owner] - 1 && index_col < mat_ranges[owner + 1]) { /* diag block */
1902           my_dnz[i] += 1;
1903         } else { /* offdiag block */
1904           my_onz[i] += 1;
1905         }
1906         /* same as before, interchanging rows and cols */
1907         if (issbaij && index_col != index_row) {
1908           owner = row_ownership[index_col];
1909           if (index_row > mat_ranges[owner] - 1 && index_row < mat_ranges[owner + 1]) {
1910             my_dnz[cols[j]] += 1;
1911           } else {
1912             my_onz[cols[j]] += 1;
1913           }
1914         }
1915       }
1916       PetscCall(MatRestoreRow(matis->A, i, &ncols, &cols, NULL));
1917     }
1918   }
1919   if (global_indices_c != global_indices_r) PetscCall(ISLocalToGlobalMappingRestoreIndices(matis->cmapping, &global_indices_c));
1920   PetscCall(ISLocalToGlobalMappingRestoreIndices(matis->rmapping, &global_indices_r));
1921   PetscCall(PetscFree(row_ownership));
1922 
1923   /* Reduce my_dnz and my_onz */
1924   if (maxreduce) {
1925     PetscCall(PetscSFReduceBegin(matis->sf, MPIU_INT, my_dnz, dnz, MPI_MAX));
1926     PetscCall(PetscSFReduceBegin(matis->sf, MPIU_INT, my_onz, onz, MPI_MAX));
1927     PetscCall(PetscSFReduceEnd(matis->sf, MPIU_INT, my_dnz, dnz, MPI_MAX));
1928     PetscCall(PetscSFReduceEnd(matis->sf, MPIU_INT, my_onz, onz, MPI_MAX));
1929   } else {
1930     PetscCall(PetscSFReduceBegin(matis->sf, MPIU_INT, my_dnz, dnz, MPI_SUM));
1931     PetscCall(PetscSFReduceBegin(matis->sf, MPIU_INT, my_onz, onz, MPI_SUM));
1932     PetscCall(PetscSFReduceEnd(matis->sf, MPIU_INT, my_dnz, dnz, MPI_SUM));
1933     PetscCall(PetscSFReduceEnd(matis->sf, MPIU_INT, my_onz, onz, MPI_SUM));
1934   }
1935   PetscCall(PetscFree2(my_dnz, my_onz));
1936 
1937   /* Resize preallocation if overestimated */
1938   for (i = 0; i < lrows; i++) {
1939     dnz[i] = PetscMin(dnz[i], lcols);
1940     onz[i] = PetscMin(onz[i], cols - lcols);
1941   }
1942 
1943   /* Set preallocation */
1944   PetscCall(MatSetBlockSizesFromMats(B, A, A));
1945   PetscCall(MatSeqAIJSetPreallocation(B, 0, dnz));
1946   PetscCall(MatMPIAIJSetPreallocation(B, 0, dnz, 0, onz));
1947   for (i = 0; i < lrows; i += bs) {
1948     PetscInt b, d = dnz[i], o = onz[i];
1949 
1950     for (b = 1; b < bs; b++) {
1951       d = PetscMax(d, dnz[i + b]);
1952       o = PetscMax(o, onz[i + b]);
1953     }
1954     dnz[i / bs] = PetscMin(d / bs + d % bs, lcols / bs);
1955     onz[i / bs] = PetscMin(o / bs + o % bs, (cols - lcols) / bs);
1956   }
1957   PetscCall(MatSeqBAIJSetPreallocation(B, bs, 0, dnz));
1958   PetscCall(MatMPIBAIJSetPreallocation(B, bs, 0, dnz, 0, onz));
1959   PetscCall(MatMPISBAIJSetPreallocation(B, bs, 0, dnz, 0, onz));
1960   MatPreallocateEnd(dnz, onz);
1961   if (issbaij) PetscCall(MatRestoreRowUpperTriangular(matis->A));
1962   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
1963   PetscFunctionReturn(PETSC_SUCCESS);
1964 }
1965 
1966 PETSC_INTERN PetscErrorCode MatConvert_IS_XAIJ(Mat mat, MatType mtype, MatReuse reuse, Mat *M)
1967 {
1968   Mat_IS            *matis = (Mat_IS *)mat->data;
1969   Mat                local_mat, MT;
1970   PetscInt           rbs, cbs, rows, cols, lrows, lcols;
1971   PetscInt           local_rows, local_cols;
1972   PetscBool          isseqdense, isseqsbaij, isseqaij, isseqbaij;
1973   PetscMPIInt        size;
1974   const PetscScalar *array;
1975 
1976   PetscFunctionBegin;
1977   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
1978   if (size == 1 && mat->rmap->N == matis->A->rmap->N && mat->cmap->N == matis->A->cmap->N && !matis->allow_repeated) {
1979     Mat      B;
1980     IS       irows = NULL, icols = NULL;
1981     PetscInt rbs, cbs;
1982 
1983     PetscCall(ISLocalToGlobalMappingGetBlockSize(matis->rmapping, &rbs));
1984     PetscCall(ISLocalToGlobalMappingGetBlockSize(matis->cmapping, &cbs));
1985     if (reuse != MAT_REUSE_MATRIX) { /* check if l2g maps are one-to-one */
1986       IS              rows, cols;
1987       const PetscInt *ridxs, *cidxs;
1988       PetscInt        i, nw;
1989       PetscBT         work;
1990 
1991       PetscCall(ISLocalToGlobalMappingGetBlockIndices(matis->rmapping, &ridxs));
1992       PetscCall(ISLocalToGlobalMappingGetSize(matis->rmapping, &nw));
1993       nw = nw / rbs;
1994       PetscCall(PetscBTCreate(nw, &work));
1995       for (i = 0; i < nw; i++) PetscCall(PetscBTSet(work, ridxs[i]));
1996       for (i = 0; i < nw; i++)
1997         if (!PetscBTLookup(work, i)) break;
1998       if (i == nw) {
1999         PetscCall(ISCreateBlock(PETSC_COMM_SELF, rbs, nw, ridxs, PETSC_USE_POINTER, &rows));
2000         PetscCall(ISSetPermutation(rows));
2001         PetscCall(ISInvertPermutation(rows, PETSC_DECIDE, &irows));
2002         PetscCall(ISDestroy(&rows));
2003       }
2004       PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(matis->rmapping, &ridxs));
2005       PetscCall(PetscBTDestroy(&work));
2006       if (irows && matis->rmapping != matis->cmapping) {
2007         PetscCall(ISLocalToGlobalMappingGetBlockIndices(matis->cmapping, &cidxs));
2008         PetscCall(ISLocalToGlobalMappingGetSize(matis->cmapping, &nw));
2009         nw = nw / cbs;
2010         PetscCall(PetscBTCreate(nw, &work));
2011         for (i = 0; i < nw; i++) PetscCall(PetscBTSet(work, cidxs[i]));
2012         for (i = 0; i < nw; i++)
2013           if (!PetscBTLookup(work, i)) break;
2014         if (i == nw) {
2015           PetscCall(ISCreateBlock(PETSC_COMM_SELF, cbs, nw, cidxs, PETSC_USE_POINTER, &cols));
2016           PetscCall(ISSetPermutation(cols));
2017           PetscCall(ISInvertPermutation(cols, PETSC_DECIDE, &icols));
2018           PetscCall(ISDestroy(&cols));
2019         }
2020         PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(matis->cmapping, &cidxs));
2021         PetscCall(PetscBTDestroy(&work));
2022       } else if (irows) {
2023         PetscCall(PetscObjectReference((PetscObject)irows));
2024         icols = irows;
2025       }
2026     } else {
2027       PetscCall(PetscObjectQuery((PetscObject)*M, "_MatIS_IS_XAIJ_irows", (PetscObject *)&irows));
2028       PetscCall(PetscObjectQuery((PetscObject)*M, "_MatIS_IS_XAIJ_icols", (PetscObject *)&icols));
2029       if (irows) PetscCall(PetscObjectReference((PetscObject)irows));
2030       if (icols) PetscCall(PetscObjectReference((PetscObject)icols));
2031     }
2032     if (!irows || !icols) {
2033       PetscCall(ISDestroy(&icols));
2034       PetscCall(ISDestroy(&irows));
2035       goto general_assembly;
2036     }
2037     PetscCall(MatConvert(matis->A, mtype, MAT_INITIAL_MATRIX, &B));
2038     if (reuse != MAT_INPLACE_MATRIX) {
2039       PetscCall(MatCreateSubMatrix(B, irows, icols, reuse, M));
2040       PetscCall(PetscObjectCompose((PetscObject)*M, "_MatIS_IS_XAIJ_irows", (PetscObject)irows));
2041       PetscCall(PetscObjectCompose((PetscObject)*M, "_MatIS_IS_XAIJ_icols", (PetscObject)icols));
2042     } else {
2043       Mat C;
2044 
2045       PetscCall(MatCreateSubMatrix(B, irows, icols, MAT_INITIAL_MATRIX, &C));
2046       PetscCall(MatHeaderReplace(mat, &C));
2047     }
2048     PetscCall(MatDestroy(&B));
2049     PetscCall(ISDestroy(&icols));
2050     PetscCall(ISDestroy(&irows));
2051     PetscFunctionReturn(PETSC_SUCCESS);
2052   }
2053 general_assembly:
2054   PetscCall(MatGetSize(mat, &rows, &cols));
2055   PetscCall(ISLocalToGlobalMappingGetBlockSize(matis->rmapping, &rbs));
2056   PetscCall(ISLocalToGlobalMappingGetBlockSize(matis->cmapping, &cbs));
2057   PetscCall(MatGetLocalSize(mat, &lrows, &lcols));
2058   PetscCall(MatGetSize(matis->A, &local_rows, &local_cols));
2059   PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A, MATSEQDENSE, &isseqdense));
2060   PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A, MATSEQAIJ, &isseqaij));
2061   PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A, MATSEQBAIJ, &isseqbaij));
2062   PetscCall(PetscObjectBaseTypeCompare((PetscObject)matis->A, MATSEQSBAIJ, &isseqsbaij));
2063   PetscCheck(isseqdense || isseqaij || isseqbaij || isseqsbaij, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not for matrix type %s", ((PetscObject)matis->A)->type_name);
2064   if (PetscDefined(USE_DEBUG)) {
2065     PetscBool lb[4], bb[4];
2066 
2067     lb[0] = isseqdense;
2068     lb[1] = isseqaij;
2069     lb[2] = isseqbaij;
2070     lb[3] = isseqsbaij;
2071     PetscCallMPI(MPIU_Allreduce(lb, bb, 4, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)mat)));
2072     PetscCheck(bb[0] || bb[1] || bb[2] || bb[3], PETSC_COMM_SELF, PETSC_ERR_SUP, "Local matrices must have the same type");
2073   }
2074 
2075   if (reuse != MAT_REUSE_MATRIX) {
2076     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &MT));
2077     PetscCall(MatSetSizes(MT, lrows, lcols, rows, cols));
2078     PetscCall(MatSetType(MT, mtype));
2079     PetscCall(MatSetBlockSizes(MT, rbs, cbs));
2080     PetscCall(MatISSetMPIXAIJPreallocation_Private(mat, MT, PETSC_FALSE));
2081   } else {
2082     PetscInt mrbs, mcbs, mrows, mcols, mlrows, mlcols;
2083 
2084     /* some checks */
2085     MT = *M;
2086     PetscCall(MatGetBlockSizes(MT, &mrbs, &mcbs));
2087     PetscCall(MatGetSize(MT, &mrows, &mcols));
2088     PetscCall(MatGetLocalSize(MT, &mlrows, &mlcols));
2089     PetscCheck(mrows == rows, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot reuse matrix. Wrong number of rows (%" PetscInt_FMT " != %" PetscInt_FMT ")", rows, mrows);
2090     PetscCheck(mcols == cols, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot reuse matrix. Wrong number of cols (%" PetscInt_FMT " != %" PetscInt_FMT ")", cols, mcols);
2091     PetscCheck(mlrows == lrows, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot reuse matrix. Wrong number of local rows (%" PetscInt_FMT " != %" PetscInt_FMT ")", lrows, mlrows);
2092     PetscCheck(mlcols == lcols, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot reuse matrix. Wrong number of local cols (%" PetscInt_FMT " != %" PetscInt_FMT ")", lcols, mlcols);
2093     PetscCheck(mrbs == rbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot reuse matrix. Wrong row block size (%" PetscInt_FMT " != %" PetscInt_FMT ")", rbs, mrbs);
2094     PetscCheck(mcbs == cbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot reuse matrix. Wrong col block size (%" PetscInt_FMT " != %" PetscInt_FMT ")", cbs, mcbs);
2095     PetscCall(MatZeroEntries(MT));
2096   }
2097 
2098   if (isseqsbaij || isseqbaij) {
2099     PetscCall(MatConvert(matis->A, MATSEQAIJ, MAT_INITIAL_MATRIX, &local_mat));
2100     isseqaij = PETSC_TRUE;
2101   } else {
2102     PetscCall(PetscObjectReference((PetscObject)matis->A));
2103     local_mat = matis->A;
2104   }
2105 
2106   /* Set values */
2107   PetscCall(MatSetLocalToGlobalMapping(MT, matis->rmapping, matis->cmapping));
2108   if (isseqdense) { /* special case for dense local matrices */
2109     PetscInt i, *dummy;
2110 
2111     PetscCall(PetscMalloc1(PetscMax(local_rows, local_cols), &dummy));
2112     for (i = 0; i < PetscMax(local_rows, local_cols); i++) dummy[i] = i;
2113     PetscCall(MatSetOption(MT, MAT_ROW_ORIENTED, PETSC_FALSE));
2114     PetscCall(MatDenseGetArrayRead(local_mat, &array));
2115     PetscCall(MatSetValuesLocal(MT, local_rows, dummy, local_cols, dummy, array, ADD_VALUES));
2116     PetscCall(MatDenseRestoreArrayRead(local_mat, &array));
2117     PetscCall(PetscFree(dummy));
2118   } else if (isseqaij) {
2119     const PetscInt *blocks;
2120     PetscInt        i, nvtxs, *xadj, *adjncy, nb;
2121     PetscBool       done;
2122     PetscScalar    *sarray;
2123 
2124     PetscCall(MatGetRowIJ(local_mat, 0, PETSC_FALSE, PETSC_FALSE, &nvtxs, (const PetscInt **)&xadj, (const PetscInt **)&adjncy, &done));
2125     PetscCheck(done, PetscObjectComm((PetscObject)local_mat), PETSC_ERR_PLIB, "Error in MatGetRowIJ");
2126     PetscCall(MatSeqAIJGetArray(local_mat, &sarray));
2127     PetscCall(MatGetVariableBlockSizes(local_mat, &nb, &blocks));
2128     if (nb) { /* speed up assembly for special blocked matrices (used by BDDC) */
2129       PetscInt sum;
2130 
2131       for (i = 0, sum = 0; i < nb; i++) sum += blocks[i];
2132       if (sum == nvtxs) {
2133         PetscInt r;
2134 
2135         for (i = 0, r = 0; i < nb; i++) {
2136           PetscAssert(blocks[i] == xadj[r + 1] - xadj[r], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid block sizes prescribed for block %" PetscInt_FMT ": expected %" PetscInt_FMT ", got %" PetscInt_FMT, i, blocks[i], xadj[r + 1] - xadj[r]);
2137           PetscCall(MatSetValuesLocal(MT, blocks[i], adjncy + xadj[r], blocks[i], adjncy + xadj[r], sarray + xadj[r], ADD_VALUES));
2138           r += blocks[i];
2139         }
2140       } else {
2141         for (i = 0; i < nvtxs; i++) PetscCall(MatSetValuesLocal(MT, 1, &i, xadj[i + 1] - xadj[i], adjncy + xadj[i], sarray + xadj[i], ADD_VALUES));
2142       }
2143     } else {
2144       for (i = 0; i < nvtxs; i++) PetscCall(MatSetValuesLocal(MT, 1, &i, xadj[i + 1] - xadj[i], PetscSafePointerPlusOffset(adjncy, xadj[i]), PetscSafePointerPlusOffset(sarray, xadj[i]), ADD_VALUES));
2145     }
2146     PetscCall(MatRestoreRowIJ(local_mat, 0, PETSC_FALSE, PETSC_FALSE, &nvtxs, (const PetscInt **)&xadj, (const PetscInt **)&adjncy, &done));
2147     PetscCheck(done, PetscObjectComm((PetscObject)local_mat), PETSC_ERR_PLIB, "Error in MatRestoreRowIJ");
2148     PetscCall(MatSeqAIJRestoreArray(local_mat, &sarray));
2149   } else { /* very basic values insertion for all other matrix types */
2150     for (PetscInt i = 0; i < local_rows; i++) {
2151       PetscInt        j;
2152       const PetscInt *local_indices_cols;
2153 
2154       PetscCall(MatGetRow(local_mat, i, &j, &local_indices_cols, &array));
2155       PetscCall(MatSetValuesLocal(MT, 1, &i, j, local_indices_cols, array, ADD_VALUES));
2156       PetscCall(MatRestoreRow(local_mat, i, &j, &local_indices_cols, &array));
2157     }
2158   }
2159   PetscCall(MatDestroy(&local_mat));
2160   PetscCall(MatAssemblyBegin(MT, MAT_FINAL_ASSEMBLY));
2161   PetscCall(MatAssemblyEnd(MT, MAT_FINAL_ASSEMBLY));
2162   if (isseqdense) PetscCall(MatSetOption(MT, MAT_ROW_ORIENTED, PETSC_TRUE));
2163   if (reuse == MAT_INPLACE_MATRIX) {
2164     PetscCall(MatHeaderReplace(mat, &MT));
2165   } else if (reuse == MAT_INITIAL_MATRIX) {
2166     *M = MT;
2167   }
2168   PetscFunctionReturn(PETSC_SUCCESS);
2169 }
2170 
2171 static PetscErrorCode MatDuplicate_IS(Mat mat, MatDuplicateOption op, Mat *newmat)
2172 {
2173   Mat_IS  *matis = (Mat_IS *)mat->data;
2174   PetscInt rbs, cbs, m, n, M, N;
2175   Mat      B, localmat;
2176 
2177   PetscFunctionBegin;
2178   PetscCall(ISLocalToGlobalMappingGetBlockSize(mat->rmap->mapping, &rbs));
2179   PetscCall(ISLocalToGlobalMappingGetBlockSize(mat->cmap->mapping, &cbs));
2180   PetscCall(MatGetSize(mat, &M, &N));
2181   PetscCall(MatGetLocalSize(mat, &m, &n));
2182   PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &B));
2183   PetscCall(MatSetSizes(B, m, n, M, N));
2184   PetscCall(MatSetBlockSize(B, rbs == cbs ? rbs : 1));
2185   PetscCall(MatSetType(B, MATIS));
2186   PetscCall(MatISSetLocalMatType(B, matis->lmattype));
2187   PetscCall(MatISSetAllowRepeated(B, matis->allow_repeated));
2188   PetscCall(MatSetLocalToGlobalMapping(B, mat->rmap->mapping, mat->cmap->mapping));
2189   PetscCall(MatDuplicate(matis->A, op, &localmat));
2190   PetscCall(MatSetLocalToGlobalMapping(localmat, matis->A->rmap->mapping, matis->A->cmap->mapping));
2191   PetscCall(MatISSetLocalMat(B, localmat));
2192   PetscCall(MatDestroy(&localmat));
2193   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2194   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2195   *newmat = B;
2196   PetscFunctionReturn(PETSC_SUCCESS);
2197 }
2198 
2199 static PetscErrorCode MatIsHermitian_IS(Mat A, PetscReal tol, PetscBool *flg)
2200 {
2201   Mat_IS   *matis = (Mat_IS *)A->data;
2202   PetscBool local_sym;
2203 
2204   PetscFunctionBegin;
2205   PetscCall(MatIsHermitian(matis->A, tol, &local_sym));
2206   PetscCallMPI(MPIU_Allreduce(&local_sym, flg, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
2207   PetscFunctionReturn(PETSC_SUCCESS);
2208 }
2209 
2210 static PetscErrorCode MatIsSymmetric_IS(Mat A, PetscReal tol, PetscBool *flg)
2211 {
2212   Mat_IS   *matis = (Mat_IS *)A->data;
2213   PetscBool local_sym;
2214 
2215   PetscFunctionBegin;
2216   if (matis->rmapping != matis->cmapping) {
2217     *flg = PETSC_FALSE;
2218     PetscFunctionReturn(PETSC_SUCCESS);
2219   }
2220   PetscCall(MatIsSymmetric(matis->A, tol, &local_sym));
2221   PetscCallMPI(MPIU_Allreduce(&local_sym, flg, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
2222   PetscFunctionReturn(PETSC_SUCCESS);
2223 }
2224 
2225 static PetscErrorCode MatIsStructurallySymmetric_IS(Mat A, PetscBool *flg)
2226 {
2227   Mat_IS   *matis = (Mat_IS *)A->data;
2228   PetscBool local_sym;
2229 
2230   PetscFunctionBegin;
2231   if (matis->rmapping != matis->cmapping) {
2232     *flg = PETSC_FALSE;
2233     PetscFunctionReturn(PETSC_SUCCESS);
2234   }
2235   PetscCall(MatIsStructurallySymmetric(matis->A, &local_sym));
2236   PetscCallMPI(MPIU_Allreduce(&local_sym, flg, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
2237   PetscFunctionReturn(PETSC_SUCCESS);
2238 }
2239 
2240 static PetscErrorCode MatDestroy_IS(Mat A)
2241 {
2242   Mat_IS *b = (Mat_IS *)A->data;
2243 
2244   PetscFunctionBegin;
2245   PetscCall(PetscFree(b->bdiag));
2246   PetscCall(PetscFree(b->lmattype));
2247   PetscCall(MatDestroy(&b->A));
2248   PetscCall(VecScatterDestroy(&b->cctx));
2249   PetscCall(VecScatterDestroy(&b->rctx));
2250   PetscCall(VecDestroy(&b->x));
2251   PetscCall(VecDestroy(&b->y));
2252   PetscCall(VecDestroy(&b->counter));
2253   PetscCall(ISDestroy(&b->getsub_ris));
2254   PetscCall(ISDestroy(&b->getsub_cis));
2255   if (b->sf != b->csf) {
2256     PetscCall(PetscSFDestroy(&b->csf));
2257     PetscCall(PetscFree2(b->csf_rootdata, b->csf_leafdata));
2258   } else b->csf = NULL;
2259   PetscCall(PetscSFDestroy(&b->sf));
2260   PetscCall(PetscFree2(b->sf_rootdata, b->sf_leafdata));
2261   PetscCall(ISLocalToGlobalMappingDestroy(&b->rmapping));
2262   PetscCall(ISLocalToGlobalMappingDestroy(&b->cmapping));
2263   PetscCall(MatDestroy(&b->dA));
2264   PetscCall(MatDestroy(&b->assembledA));
2265   PetscCall(PetscFree(A->data));
2266   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
2267   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISSetLocalMatType_C", NULL));
2268   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISGetLocalMat_C", NULL));
2269   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISSetLocalMat_C", NULL));
2270   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISRestoreLocalMat_C", NULL));
2271   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISSetPreallocation_C", NULL));
2272   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISStoreL2L_C", NULL));
2273   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISFixLocalEmpty_C", NULL));
2274   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISGetLocalToGlobalMapping_C", NULL));
2275   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_mpiaij_C", NULL));
2276   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_mpibaij_C", NULL));
2277   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_mpisbaij_C", NULL));
2278   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_seqaij_C", NULL));
2279   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_seqbaij_C", NULL));
2280   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_seqsbaij_C", NULL));
2281   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_aij_C", NULL));
2282   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", NULL));
2283   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetPreallocationCOO_C", NULL));
2284   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetValuesCOO_C", NULL));
2285   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISSetAllowRepeated_C", NULL));
2286   PetscFunctionReturn(PETSC_SUCCESS);
2287 }
2288 
2289 static PetscErrorCode MatMult_IS(Mat A, Vec x, Vec y)
2290 {
2291   Mat_IS     *is   = (Mat_IS *)A->data;
2292   PetscScalar zero = 0.0;
2293 
2294   PetscFunctionBegin;
2295   /*  scatter the global vector x into the local work vector */
2296   PetscCall(VecScatterBegin(is->cctx, x, is->x, INSERT_VALUES, SCATTER_FORWARD));
2297   PetscCall(VecScatterEnd(is->cctx, x, is->x, INSERT_VALUES, SCATTER_FORWARD));
2298 
2299   /* multiply the local matrix */
2300   PetscCall(MatMult(is->A, is->x, is->y));
2301 
2302   /* scatter product back into global memory */
2303   PetscCall(VecSet(y, zero));
2304   PetscCall(VecScatterBegin(is->rctx, is->y, y, ADD_VALUES, SCATTER_REVERSE));
2305   PetscCall(VecScatterEnd(is->rctx, is->y, y, ADD_VALUES, SCATTER_REVERSE));
2306   PetscFunctionReturn(PETSC_SUCCESS);
2307 }
2308 
2309 static PetscErrorCode MatMultAdd_IS(Mat A, Vec v1, Vec v2, Vec v3)
2310 {
2311   Vec temp_vec;
2312 
2313   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
2314   if (v3 != v2) {
2315     PetscCall(MatMult(A, v1, v3));
2316     PetscCall(VecAXPY(v3, 1.0, v2));
2317   } else {
2318     PetscCall(VecDuplicate(v2, &temp_vec));
2319     PetscCall(MatMult(A, v1, temp_vec));
2320     PetscCall(VecAXPY(temp_vec, 1.0, v2));
2321     PetscCall(VecCopy(temp_vec, v3));
2322     PetscCall(VecDestroy(&temp_vec));
2323   }
2324   PetscFunctionReturn(PETSC_SUCCESS);
2325 }
2326 
2327 static PetscErrorCode MatMultTranspose_IS(Mat A, Vec y, Vec x)
2328 {
2329   Mat_IS *is = (Mat_IS *)A->data;
2330 
2331   PetscFunctionBegin;
2332   /*  scatter the global vector x into the local work vector */
2333   PetscCall(VecScatterBegin(is->rctx, y, is->y, INSERT_VALUES, SCATTER_FORWARD));
2334   PetscCall(VecScatterEnd(is->rctx, y, is->y, INSERT_VALUES, SCATTER_FORWARD));
2335 
2336   /* multiply the local matrix */
2337   PetscCall(MatMultTranspose(is->A, is->y, is->x));
2338 
2339   /* scatter product back into global vector */
2340   PetscCall(VecSet(x, 0));
2341   PetscCall(VecScatterBegin(is->cctx, is->x, x, ADD_VALUES, SCATTER_REVERSE));
2342   PetscCall(VecScatterEnd(is->cctx, is->x, x, ADD_VALUES, SCATTER_REVERSE));
2343   PetscFunctionReturn(PETSC_SUCCESS);
2344 }
2345 
2346 static PetscErrorCode MatMultTransposeAdd_IS(Mat A, Vec v1, Vec v2, Vec v3)
2347 {
2348   Vec temp_vec;
2349 
2350   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
2351   if (v3 != v2) {
2352     PetscCall(MatMultTranspose(A, v1, v3));
2353     PetscCall(VecAXPY(v3, 1.0, v2));
2354   } else {
2355     PetscCall(VecDuplicate(v2, &temp_vec));
2356     PetscCall(MatMultTranspose(A, v1, temp_vec));
2357     PetscCall(VecAXPY(temp_vec, 1.0, v2));
2358     PetscCall(VecCopy(temp_vec, v3));
2359     PetscCall(VecDestroy(&temp_vec));
2360   }
2361   PetscFunctionReturn(PETSC_SUCCESS);
2362 }
2363 
2364 static PetscErrorCode MatView_IS(Mat A, PetscViewer viewer)
2365 {
2366   Mat_IS                *a = (Mat_IS *)A->data;
2367   PetscViewer            sviewer;
2368   PetscBool              isascii, isbinary, viewl2g = PETSC_FALSE, native;
2369   PetscViewerFormat      format;
2370   ISLocalToGlobalMapping rmap, cmap;
2371 
2372   PetscFunctionBegin;
2373   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2374   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
2375   PetscCall(PetscViewerGetFormat(viewer, &format));
2376   native = (PetscBool)(format == PETSC_VIEWER_NATIVE);
2377   if (native) {
2378     rmap = A->rmap->mapping;
2379     cmap = A->cmap->mapping;
2380   } else {
2381     rmap = a->rmapping;
2382     cmap = a->cmapping;
2383   }
2384   if (isascii) {
2385     if (format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS);
2386     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_MATLAB) viewl2g = PETSC_TRUE;
2387   } else if (isbinary) {
2388     PetscInt    tr[6], nr, nc;
2389     char        lmattype[64] = {'\0'};
2390     PetscMPIInt size;
2391     PetscBool   skipHeader;
2392     IS          is;
2393 
2394     PetscCall(PetscViewerSetUp(viewer));
2395     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
2396     tr[0] = MAT_FILE_CLASSID;
2397     tr[1] = A->rmap->N;
2398     tr[2] = A->cmap->N;
2399     tr[3] = -size; /* AIJ stores nnz here */
2400     tr[4] = (PetscInt)(rmap == cmap);
2401     tr[5] = a->allow_repeated;
2402     PetscCall(PetscSNPrintf(lmattype, sizeof(lmattype), "%s", a->lmattype));
2403 
2404     PetscCall(PetscViewerBinaryWrite(viewer, tr, PETSC_STATIC_ARRAY_LENGTH(tr), PETSC_INT));
2405     PetscCall(PetscViewerBinaryWrite(viewer, lmattype, sizeof(lmattype), PETSC_CHAR));
2406 
2407     /* first dump l2g info (we need the header for proper loading on different number of processes) */
2408     PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
2409     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_FALSE));
2410     PetscCall(ISLocalToGlobalMappingView(rmap, viewer));
2411     if (cmap != rmap) PetscCall(ISLocalToGlobalMappingView(cmap, viewer));
2412 
2413     /* then the sizes of the local matrices */
2414     PetscCall(MatGetSize(a->A, &nr, &nc));
2415     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)viewer), 1, &nr, PETSC_USE_POINTER, &is));
2416     PetscCall(ISView(is, viewer));
2417     PetscCall(ISDestroy(&is));
2418     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)viewer), 1, &nc, PETSC_USE_POINTER, &is));
2419     PetscCall(ISView(is, viewer));
2420     PetscCall(ISDestroy(&is));
2421     PetscCall(PetscViewerBinarySetSkipHeader(viewer, skipHeader));
2422   }
2423   if (format == PETSC_VIEWER_ASCII_MATLAB) {
2424     char        name[64];
2425     PetscMPIInt size, rank;
2426 
2427     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
2428     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
2429     if (size > 1) PetscCall(PetscSNPrintf(name, sizeof(name), "lmat_%d", rank));
2430     else PetscCall(PetscSNPrintf(name, sizeof(name), "lmat"));
2431     PetscCall(PetscObjectSetName((PetscObject)a->A, name));
2432   }
2433 
2434   /* Dump the local matrices */
2435   if (isbinary) { /* ViewerGetSubViewer does not work in parallel */
2436     PetscBool   isaij;
2437     PetscInt    nr, nc;
2438     Mat         lA, B;
2439     Mat_MPIAIJ *b;
2440 
2441     /* We create a temporary MPIAIJ matrix that stores the unassembled operator */
2442     PetscCall(PetscObjectBaseTypeCompare((PetscObject)a->A, MATAIJ, &isaij));
2443     if (!isaij) PetscCall(MatConvert(a->A, MATSEQAIJ, MAT_INITIAL_MATRIX, &lA));
2444     else {
2445       PetscCall(PetscObjectReference((PetscObject)a->A));
2446       lA = a->A;
2447     }
2448     PetscCall(MatCreate(PetscObjectComm((PetscObject)viewer), &B));
2449     PetscCall(MatSetType(B, MATMPIAIJ));
2450     PetscCall(MatGetSize(lA, &nr, &nc));
2451     PetscCall(MatSetSizes(B, nr, nc, PETSC_DECIDE, PETSC_DECIDE));
2452     PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL));
2453 
2454     b = (Mat_MPIAIJ *)B->data;
2455     PetscCall(MatDestroy(&b->A));
2456     b->A = lA;
2457 
2458     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2459     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2460     PetscCall(MatView(B, viewer));
2461     PetscCall(MatDestroy(&B));
2462   } else {
2463     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2464     PetscCall(MatView(a->A, sviewer));
2465     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2466   }
2467 
2468   /* with ASCII, we dump the l2gmaps at the end */
2469   if (viewl2g) {
2470     if (format == PETSC_VIEWER_ASCII_MATLAB) {
2471       PetscCall(PetscObjectSetName((PetscObject)rmap, "row"));
2472       PetscCall(ISLocalToGlobalMappingView(rmap, viewer));
2473       PetscCall(PetscObjectSetName((PetscObject)cmap, "col"));
2474       PetscCall(ISLocalToGlobalMappingView(cmap, viewer));
2475     } else {
2476       PetscCall(ISLocalToGlobalMappingView(rmap, viewer));
2477       PetscCall(ISLocalToGlobalMappingView(cmap, viewer));
2478     }
2479   }
2480   PetscFunctionReturn(PETSC_SUCCESS);
2481 }
2482 
2483 static PetscErrorCode MatLoad_IS(Mat A, PetscViewer viewer)
2484 {
2485   ISLocalToGlobalMapping rmap, cmap;
2486   MPI_Comm               comm = PetscObjectComm((PetscObject)A);
2487   PetscBool              isbinary, samel, allow, isbaij;
2488   PetscInt               tr[6], M, N, nr, nc, Asize, isn;
2489   const PetscInt        *idx;
2490   PetscMPIInt            size;
2491   char                   lmattype[64];
2492   Mat                    dA, lA;
2493   IS                     is;
2494 
2495   PetscFunctionBegin;
2496   PetscCheckSameComm(A, 1, viewer, 2);
2497   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
2498   PetscCheck(isbinary, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Invalid viewer of type %s", ((PetscObject)viewer)->type_name);
2499 
2500   PetscCall(PetscViewerBinaryRead(viewer, tr, PETSC_STATIC_ARRAY_LENGTH(tr), NULL, PETSC_INT));
2501   PetscCheck(tr[0] == MAT_FILE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a matrix next in file");
2502   PetscCheck(tr[1] >= 0, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a IS matrix next in file");
2503   PetscCheck(tr[2] >= 0, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a IS matrix next in file");
2504   PetscCheck(tr[3] < 0, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a IS matrix next in file");
2505   PetscCheck(tr[4] == 0 || tr[4] == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a IS matrix next in file");
2506   PetscCheck(tr[5] == 0 || tr[5] == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a IS matrix next in file");
2507   M     = tr[1];
2508   N     = tr[2];
2509   Asize = -tr[3];
2510   samel = (PetscBool)tr[4];
2511   allow = (PetscBool)tr[5];
2512   PetscCall(PetscViewerBinaryRead(viewer, lmattype, sizeof(lmattype), NULL, PETSC_CHAR));
2513 
2514   /* if we are loading from a larger set of processes, allow repeated entries */
2515   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
2516   if (Asize > size) allow = PETSC_TRUE;
2517 
2518   /* set global sizes if not set already */
2519   if (A->rmap->N < 0) A->rmap->N = M;
2520   if (A->cmap->N < 0) A->cmap->N = N;
2521   PetscCall(PetscLayoutSetUp(A->rmap));
2522   PetscCall(PetscLayoutSetUp(A->cmap));
2523   PetscCheck(M == A->rmap->N, comm, PETSC_ERR_ARG_SIZ, "Matrix rows should be %" PetscInt_FMT ", found %" PetscInt_FMT, M, A->rmap->N);
2524   PetscCheck(N == A->cmap->N, comm, PETSC_ERR_ARG_SIZ, "Matrix columns should be %" PetscInt_FMT ", found %" PetscInt_FMT, N, A->cmap->N);
2525 
2526   /* load l2g maps */
2527   PetscCall(ISLocalToGlobalMappingCreate(comm, 0, 0, NULL, PETSC_USE_POINTER, &rmap));
2528   PetscCall(ISLocalToGlobalMappingLoad(rmap, viewer));
2529   if (!samel) {
2530     PetscCall(ISLocalToGlobalMappingCreate(comm, 0, 0, NULL, PETSC_USE_POINTER, &cmap));
2531     PetscCall(ISLocalToGlobalMappingLoad(cmap, viewer));
2532   } else {
2533     PetscCall(PetscObjectReference((PetscObject)rmap));
2534     cmap = rmap;
2535   }
2536 
2537   /* load sizes of local matrices */
2538   PetscCall(ISCreate(comm, &is));
2539   PetscCall(ISSetType(is, ISGENERAL));
2540   PetscCall(ISLoad(is, viewer));
2541   PetscCall(ISGetLocalSize(is, &isn));
2542   PetscCall(ISGetIndices(is, &idx));
2543   nr = 0;
2544   for (PetscInt i = 0; i < isn; i++) nr += idx[i];
2545   PetscCall(ISRestoreIndices(is, &idx));
2546   PetscCall(ISDestroy(&is));
2547   PetscCall(ISCreate(comm, &is));
2548   PetscCall(ISSetType(is, ISGENERAL));
2549   PetscCall(ISLoad(is, viewer));
2550   PetscCall(ISGetLocalSize(is, &isn));
2551   PetscCall(ISGetIndices(is, &idx));
2552   nc = 0;
2553   for (PetscInt i = 0; i < isn; i++) nc += idx[i];
2554   PetscCall(ISRestoreIndices(is, &idx));
2555   PetscCall(ISDestroy(&is));
2556 
2557   /* now load the unassembled operator */
2558   PetscCall(MatCreate(comm, &dA));
2559   PetscCall(MatSetType(dA, MATMPIAIJ));
2560   PetscCall(MatSetSizes(dA, nr, nc, PETSC_DECIDE, PETSC_DECIDE));
2561   PetscCall(MatLoad(dA, viewer));
2562   PetscCall(MatMPIAIJGetSeqAIJ(dA, &lA, NULL, NULL));
2563   PetscCall(PetscObjectReference((PetscObject)lA));
2564   PetscCall(MatDestroy(&dA));
2565 
2566   /* and convert to the desired format */
2567   PetscCall(PetscStrcmpAny(lmattype, &isbaij, MATSBAIJ, MATSEQSBAIJ, ""));
2568   if (isbaij) PetscCall(MatSetOption(lA, MAT_SYMMETRIC, PETSC_TRUE));
2569   PetscCall(MatConvert(lA, lmattype, MAT_INPLACE_MATRIX, &lA));
2570 
2571   /* assemble the MATIS object */
2572   PetscCall(MatISSetAllowRepeated(A, allow));
2573   PetscCall(MatSetLocalToGlobalMapping(A, rmap, cmap));
2574   PetscCall(MatISSetLocalMat(A, lA));
2575   PetscCall(MatDestroy(&lA));
2576   PetscCall(ISLocalToGlobalMappingDestroy(&rmap));
2577   PetscCall(ISLocalToGlobalMappingDestroy(&cmap));
2578   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2579   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2580   PetscFunctionReturn(PETSC_SUCCESS);
2581 }
2582 
2583 static PetscErrorCode MatInvertBlockDiagonal_IS(Mat mat, const PetscScalar **values)
2584 {
2585   Mat_IS            *is = (Mat_IS *)mat->data;
2586   MPI_Datatype       nodeType;
2587   const PetscScalar *lv;
2588   PetscInt           bs;
2589   PetscMPIInt        mbs;
2590 
2591   PetscFunctionBegin;
2592   PetscCall(MatGetBlockSize(mat, &bs));
2593   PetscCall(MatSetBlockSize(is->A, bs));
2594   PetscCall(MatInvertBlockDiagonal(is->A, &lv));
2595   if (!is->bdiag) PetscCall(PetscMalloc1(bs * mat->rmap->n, &is->bdiag));
2596   PetscCall(PetscMPIIntCast(bs, &mbs));
2597   PetscCallMPI(MPI_Type_contiguous(mbs, MPIU_SCALAR, &nodeType));
2598   PetscCallMPI(MPI_Type_commit(&nodeType));
2599   PetscCall(PetscSFReduceBegin(is->sf, nodeType, lv, is->bdiag, MPI_REPLACE));
2600   PetscCall(PetscSFReduceEnd(is->sf, nodeType, lv, is->bdiag, MPI_REPLACE));
2601   PetscCallMPI(MPI_Type_free(&nodeType));
2602   if (values) *values = is->bdiag;
2603   PetscFunctionReturn(PETSC_SUCCESS);
2604 }
2605 
2606 static PetscErrorCode MatISSetUpScatters_Private(Mat A)
2607 {
2608   Vec             cglobal, rglobal;
2609   IS              from;
2610   Mat_IS         *is = (Mat_IS *)A->data;
2611   PetscScalar     sum;
2612   const PetscInt *garray;
2613   PetscInt        nr, rbs, nc, cbs;
2614   VecType         rtype;
2615 
2616   PetscFunctionBegin;
2617   PetscCall(ISLocalToGlobalMappingGetSize(is->rmapping, &nr));
2618   PetscCall(ISLocalToGlobalMappingGetBlockSize(is->rmapping, &rbs));
2619   PetscCall(ISLocalToGlobalMappingGetSize(is->cmapping, &nc));
2620   PetscCall(ISLocalToGlobalMappingGetBlockSize(is->cmapping, &cbs));
2621   PetscCall(VecDestroy(&is->x));
2622   PetscCall(VecDestroy(&is->y));
2623   PetscCall(VecDestroy(&is->counter));
2624   PetscCall(VecScatterDestroy(&is->rctx));
2625   PetscCall(VecScatterDestroy(&is->cctx));
2626   PetscCall(MatCreateVecs(is->A, &is->x, &is->y));
2627   PetscCall(VecBindToCPU(is->y, PETSC_TRUE));
2628   PetscCall(VecGetRootType_Private(is->y, &rtype));
2629   PetscCall(PetscFree(A->defaultvectype));
2630   PetscCall(PetscStrallocpy(rtype, &A->defaultvectype));
2631   PetscCall(MatCreateVecs(A, &cglobal, &rglobal));
2632   PetscCall(VecBindToCPU(rglobal, PETSC_TRUE));
2633   PetscCall(ISLocalToGlobalMappingGetBlockIndices(is->rmapping, &garray));
2634   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A), rbs, nr / rbs, garray, PETSC_USE_POINTER, &from));
2635   PetscCall(VecScatterCreate(rglobal, from, is->y, NULL, &is->rctx));
2636   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(is->rmapping, &garray));
2637   PetscCall(ISDestroy(&from));
2638   if (is->rmapping != is->cmapping) {
2639     PetscCall(ISLocalToGlobalMappingGetBlockIndices(is->cmapping, &garray));
2640     PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A), cbs, nc / cbs, garray, PETSC_USE_POINTER, &from));
2641     PetscCall(VecScatterCreate(cglobal, from, is->x, NULL, &is->cctx));
2642     PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(is->cmapping, &garray));
2643     PetscCall(ISDestroy(&from));
2644   } else {
2645     PetscCall(PetscObjectReference((PetscObject)is->rctx));
2646     is->cctx = is->rctx;
2647   }
2648   PetscCall(VecDestroy(&cglobal));
2649 
2650   /* interface counter vector (local) */
2651   PetscCall(VecDuplicate(is->y, &is->counter));
2652   PetscCall(VecBindToCPU(is->counter, PETSC_TRUE));
2653   PetscCall(VecSet(is->y, 1.));
2654   PetscCall(VecScatterBegin(is->rctx, is->y, rglobal, ADD_VALUES, SCATTER_REVERSE));
2655   PetscCall(VecScatterEnd(is->rctx, is->y, rglobal, ADD_VALUES, SCATTER_REVERSE));
2656   PetscCall(VecScatterBegin(is->rctx, rglobal, is->counter, INSERT_VALUES, SCATTER_FORWARD));
2657   PetscCall(VecScatterEnd(is->rctx, rglobal, is->counter, INSERT_VALUES, SCATTER_FORWARD));
2658   PetscCall(VecBindToCPU(is->y, PETSC_FALSE));
2659   PetscCall(VecBindToCPU(is->counter, PETSC_FALSE));
2660 
2661   /* special functions for block-diagonal matrices */
2662   PetscCall(VecSum(rglobal, &sum));
2663   A->ops->invertblockdiagonal = NULL;
2664   if ((PetscInt)(PetscRealPart(sum)) == A->rmap->N && A->rmap->N == A->cmap->N && is->rmapping == is->cmapping) A->ops->invertblockdiagonal = MatInvertBlockDiagonal_IS;
2665   PetscCall(VecDestroy(&rglobal));
2666 
2667   /* setup SF for general purpose shared indices based communications */
2668   PetscCall(MatISSetUpSF_IS(A));
2669   PetscFunctionReturn(PETSC_SUCCESS);
2670 }
2671 
2672 static PetscErrorCode MatISFilterL2GMap(Mat A, ISLocalToGlobalMapping map, ISLocalToGlobalMapping *nmap, ISLocalToGlobalMapping *lmap)
2673 {
2674   Mat_IS                    *matis = (Mat_IS *)A->data;
2675   IS                         is;
2676   ISLocalToGlobalMappingType l2gtype;
2677   const PetscInt            *idxs;
2678   PetscHSetI                 ht;
2679   PetscInt                  *nidxs;
2680   PetscInt                   i, n, bs, c;
2681   PetscBool                  flg[] = {PETSC_FALSE, PETSC_FALSE};
2682 
2683   PetscFunctionBegin;
2684   PetscCall(ISLocalToGlobalMappingGetSize(map, &n));
2685   PetscCall(ISLocalToGlobalMappingGetBlockSize(map, &bs));
2686   PetscCall(ISLocalToGlobalMappingGetBlockIndices(map, &idxs));
2687   PetscCall(PetscHSetICreate(&ht));
2688   PetscCall(PetscMalloc1(n / bs, &nidxs));
2689   for (i = 0, c = 0; i < n / bs; i++) {
2690     PetscBool missing = PETSC_TRUE;
2691     if (idxs[i] < 0) {
2692       flg[0] = PETSC_TRUE;
2693       continue;
2694     }
2695     if (!matis->allow_repeated) PetscCall(PetscHSetIQueryAdd(ht, idxs[i], &missing));
2696     if (!missing) flg[1] = PETSC_TRUE;
2697     else nidxs[c++] = idxs[i];
2698   }
2699   PetscCall(PetscHSetIDestroy(&ht));
2700   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, flg, 2, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)A)));
2701   if (!flg[0] && !flg[1]) { /* Entries are all non negative and unique */
2702     *nmap = NULL;
2703     *lmap = NULL;
2704     PetscCall(PetscFree(nidxs));
2705     PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(map, &idxs));
2706     PetscFunctionReturn(PETSC_SUCCESS);
2707   }
2708 
2709   /* New l2g map without negative indices (and repeated indices if not allowed) */
2710   PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A), bs, c, nidxs, PETSC_USE_POINTER, &is));
2711   PetscCall(ISLocalToGlobalMappingCreateIS(is, nmap));
2712   PetscCall(ISDestroy(&is));
2713   PetscCall(ISLocalToGlobalMappingGetType(map, &l2gtype));
2714   PetscCall(ISLocalToGlobalMappingSetType(*nmap, l2gtype));
2715 
2716   /* New local l2g map for repeated indices if not allowed */
2717   PetscCall(ISGlobalToLocalMappingApplyBlock(*nmap, IS_GTOLM_MASK, n / bs, idxs, NULL, nidxs));
2718   PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, n / bs, nidxs, PETSC_USE_POINTER, &is));
2719   PetscCall(ISLocalToGlobalMappingCreateIS(is, lmap));
2720   PetscCall(ISDestroy(&is));
2721   PetscCall(PetscFree(nidxs));
2722   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(map, &idxs));
2723   PetscFunctionReturn(PETSC_SUCCESS);
2724 }
2725 
2726 static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A, ISLocalToGlobalMapping rmapping, ISLocalToGlobalMapping cmapping)
2727 {
2728   Mat_IS                *is            = (Mat_IS *)A->data;
2729   ISLocalToGlobalMapping localrmapping = NULL, localcmapping = NULL;
2730   PetscInt               nr, rbs, nc, cbs;
2731   PetscBool              cong, freem[] = {PETSC_FALSE, PETSC_FALSE};
2732 
2733   PetscFunctionBegin;
2734   if (rmapping) PetscCheckSameComm(A, 1, rmapping, 2);
2735   if (cmapping) PetscCheckSameComm(A, 1, cmapping, 3);
2736 
2737   PetscCall(ISLocalToGlobalMappingDestroy(&is->rmapping));
2738   PetscCall(ISLocalToGlobalMappingDestroy(&is->cmapping));
2739   PetscCall(PetscLayoutSetUp(A->rmap));
2740   PetscCall(PetscLayoutSetUp(A->cmap));
2741   PetscCall(MatHasCongruentLayouts(A, &cong));
2742 
2743   /* If NULL, local space matches global space */
2744   if (!rmapping) {
2745     IS is;
2746 
2747     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A), A->rmap->N, 0, 1, &is));
2748     PetscCall(ISLocalToGlobalMappingCreateIS(is, &rmapping));
2749     PetscCall(ISLocalToGlobalMappingSetBlockSize(rmapping, A->rmap->bs));
2750     PetscCall(ISDestroy(&is));
2751     freem[0] = PETSC_TRUE;
2752     if (!cmapping && cong && A->rmap->bs == A->cmap->bs) cmapping = rmapping;
2753   } else if (!is->islocalref) { /* check if the l2g map has negative or repeated entries */
2754     PetscCall(MatISFilterL2GMap(A, rmapping, &is->rmapping, &localrmapping));
2755     if (rmapping == cmapping) {
2756       PetscCall(PetscObjectReference((PetscObject)is->rmapping));
2757       is->cmapping = is->rmapping;
2758       PetscCall(PetscObjectReference((PetscObject)localrmapping));
2759       localcmapping = localrmapping;
2760     }
2761   }
2762   if (!cmapping) {
2763     IS is;
2764 
2765     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A), A->cmap->N, 0, 1, &is));
2766     PetscCall(ISLocalToGlobalMappingCreateIS(is, &cmapping));
2767     PetscCall(ISLocalToGlobalMappingSetBlockSize(cmapping, A->cmap->bs));
2768     PetscCall(ISDestroy(&is));
2769     freem[1] = PETSC_TRUE;
2770   } else if (cmapping != rmapping && !is->islocalref) { /* check if the l2g map has negative or repeated entries */
2771     PetscCall(MatISFilterL2GMap(A, cmapping, &is->cmapping, &localcmapping));
2772   }
2773   if (!is->rmapping) {
2774     PetscCall(PetscObjectReference((PetscObject)rmapping));
2775     is->rmapping = rmapping;
2776   }
2777   if (!is->cmapping) {
2778     PetscCall(PetscObjectReference((PetscObject)cmapping));
2779     is->cmapping = cmapping;
2780   }
2781 
2782   /* Clean up */
2783   is->lnnzstate = 0;
2784   PetscCall(MatDestroy(&is->dA));
2785   PetscCall(MatDestroy(&is->assembledA));
2786   PetscCall(MatDestroy(&is->A));
2787   if (is->csf != is->sf) {
2788     PetscCall(PetscSFDestroy(&is->csf));
2789     PetscCall(PetscFree2(is->csf_rootdata, is->csf_leafdata));
2790   } else is->csf = NULL;
2791   PetscCall(PetscSFDestroy(&is->sf));
2792   PetscCall(PetscFree2(is->sf_rootdata, is->sf_leafdata));
2793   PetscCall(PetscFree(is->bdiag));
2794 
2795   /* check if the two mappings are actually the same for square matrices since MATIS has some optimization for this case
2796      (DOLFIN passes 2 different objects) */
2797   PetscCall(ISLocalToGlobalMappingGetSize(is->rmapping, &nr));
2798   PetscCall(ISLocalToGlobalMappingGetBlockSize(is->rmapping, &rbs));
2799   PetscCall(ISLocalToGlobalMappingGetSize(is->cmapping, &nc));
2800   PetscCall(ISLocalToGlobalMappingGetBlockSize(is->cmapping, &cbs));
2801   if (is->rmapping != is->cmapping && cong) {
2802     PetscBool same = PETSC_FALSE;
2803     if (nr == nc && cbs == rbs) {
2804       const PetscInt *idxs1, *idxs2;
2805 
2806       PetscCall(ISLocalToGlobalMappingGetBlockIndices(is->rmapping, &idxs1));
2807       PetscCall(ISLocalToGlobalMappingGetBlockIndices(is->cmapping, &idxs2));
2808       PetscCall(PetscArraycmp(idxs1, idxs2, nr / rbs, &same));
2809       PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(is->rmapping, &idxs1));
2810       PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(is->cmapping, &idxs2));
2811     }
2812     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &same, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
2813     if (same) {
2814       PetscCall(ISLocalToGlobalMappingDestroy(&is->cmapping));
2815       PetscCall(PetscObjectReference((PetscObject)is->rmapping));
2816       is->cmapping = is->rmapping;
2817     }
2818   }
2819   PetscCall(PetscLayoutSetBlockSize(A->rmap, rbs));
2820   PetscCall(PetscLayoutSetBlockSize(A->cmap, cbs));
2821   /* Pass the user defined maps to the layout */
2822   PetscCall(PetscLayoutSetISLocalToGlobalMapping(A->rmap, rmapping));
2823   PetscCall(PetscLayoutSetISLocalToGlobalMapping(A->cmap, cmapping));
2824   if (freem[0]) PetscCall(ISLocalToGlobalMappingDestroy(&rmapping));
2825   if (freem[1]) PetscCall(ISLocalToGlobalMappingDestroy(&cmapping));
2826 
2827   /* Create the local matrix A */
2828   PetscCall(MatCreate(PETSC_COMM_SELF, &is->A));
2829   PetscCall(MatSetType(is->A, is->lmattype));
2830   PetscCall(MatSetSizes(is->A, nr, nc, nr, nc));
2831   PetscCall(MatSetBlockSizes(is->A, rbs, cbs));
2832   PetscCall(MatSetOptionsPrefix(is->A, "is_"));
2833   PetscCall(MatAppendOptionsPrefix(is->A, ((PetscObject)A)->prefix));
2834   PetscCall(PetscLayoutSetUp(is->A->rmap));
2835   PetscCall(PetscLayoutSetUp(is->A->cmap));
2836   PetscCall(MatSetLocalToGlobalMapping(is->A, localrmapping, localcmapping));
2837   PetscCall(ISLocalToGlobalMappingDestroy(&localrmapping));
2838   PetscCall(ISLocalToGlobalMappingDestroy(&localcmapping));
2839 
2840   /* setup scatters and local vectors for MatMult */
2841   if (!is->islocalref) PetscCall(MatISSetUpScatters_Private(A));
2842   A->preallocated = PETSC_TRUE;
2843   PetscFunctionReturn(PETSC_SUCCESS);
2844 }
2845 
2846 static PetscErrorCode MatSetUp_IS(Mat A)
2847 {
2848   Mat_IS                *is = (Mat_IS *)A->data;
2849   ISLocalToGlobalMapping rmap, cmap;
2850 
2851   PetscFunctionBegin;
2852   if (!is->sf) {
2853     PetscCall(MatGetLocalToGlobalMapping(A, &rmap, &cmap));
2854     PetscCall(MatSetLocalToGlobalMapping(A, rmap, cmap));
2855   }
2856   PetscFunctionReturn(PETSC_SUCCESS);
2857 }
2858 
2859 static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode addv)
2860 {
2861   Mat_IS  *is = (Mat_IS *)mat->data;
2862   PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION], cols_l[MATIS_MAX_ENTRIES_INSERTION];
2863 
2864   PetscFunctionBegin;
2865   PetscCall(ISGlobalToLocalMappingApply(is->rmapping, IS_GTOLM_MASK, m, rows, &m, rows_l));
2866   if (m != n || rows != cols || is->cmapping != is->rmapping) {
2867     PetscCall(ISGlobalToLocalMappingApply(is->cmapping, IS_GTOLM_MASK, n, cols, &n, cols_l));
2868     PetscCall(MatSetValues(is->A, m, rows_l, n, cols_l, values, addv));
2869   } else {
2870     PetscCall(MatSetValues(is->A, m, rows_l, m, rows_l, values, addv));
2871   }
2872   PetscFunctionReturn(PETSC_SUCCESS);
2873 }
2874 
2875 static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode addv)
2876 {
2877   Mat_IS  *is = (Mat_IS *)mat->data;
2878   PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION], cols_l[MATIS_MAX_ENTRIES_INSERTION];
2879 
2880   PetscFunctionBegin;
2881   PetscCall(ISGlobalToLocalMappingApplyBlock(is->rmapping, IS_GTOLM_MASK, m, rows, &m, rows_l));
2882   if (m != n || rows != cols || is->cmapping != is->rmapping) {
2883     PetscCall(ISGlobalToLocalMappingApplyBlock(is->cmapping, IS_GTOLM_MASK, n, cols, &n, cols_l));
2884     PetscCall(MatSetValuesBlocked(is->A, m, rows_l, n, cols_l, values, addv));
2885   } else {
2886     PetscCall(MatSetValuesBlocked(is->A, m, rows_l, m, rows_l, values, addv));
2887   }
2888   PetscFunctionReturn(PETSC_SUCCESS);
2889 }
2890 
2891 static PetscErrorCode MatSetValuesLocal_IS(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode addv)
2892 {
2893   Mat_IS *is = (Mat_IS *)A->data;
2894 
2895   PetscFunctionBegin;
2896   if (is->A->rmap->mapping || is->A->cmap->mapping) {
2897     PetscCall(MatSetValuesLocal(is->A, m, rows, n, cols, values, addv));
2898   } else {
2899     PetscCall(MatSetValues(is->A, m, rows, n, cols, values, addv));
2900   }
2901   PetscFunctionReturn(PETSC_SUCCESS);
2902 }
2903 
2904 static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode addv)
2905 {
2906   Mat_IS *is = (Mat_IS *)A->data;
2907 
2908   PetscFunctionBegin;
2909   if (is->A->rmap->mapping || is->A->cmap->mapping) {
2910     PetscCall(MatSetValuesBlockedLocal(is->A, m, rows, n, cols, values, addv));
2911   } else {
2912     PetscCall(MatSetValuesBlocked(is->A, m, rows, n, cols, values, addv));
2913   }
2914   PetscFunctionReturn(PETSC_SUCCESS);
2915 }
2916 
2917 static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A, PetscInt n, const PetscInt rows[], PetscScalar diag, PetscBool columns)
2918 {
2919   Mat_IS *is = (Mat_IS *)A->data;
2920 
2921   PetscFunctionBegin;
2922   if (!n) {
2923     is->pure_neumann = PETSC_TRUE;
2924   } else {
2925     PetscInt i;
2926     is->pure_neumann = PETSC_FALSE;
2927 
2928     if (columns) {
2929       PetscCall(MatZeroRowsColumns(is->A, n, rows, diag, NULL, NULL));
2930     } else {
2931       PetscCall(MatZeroRows(is->A, n, rows, diag, NULL, NULL));
2932     }
2933     if (diag != 0.) {
2934       const PetscScalar *array;
2935       PetscCall(VecGetArrayRead(is->counter, &array));
2936       for (i = 0; i < n; i++) PetscCall(MatSetValue(is->A, rows[i], rows[i], diag / (array[rows[i]]), INSERT_VALUES));
2937       PetscCall(VecRestoreArrayRead(is->counter, &array));
2938     }
2939     PetscCall(MatAssemblyBegin(is->A, MAT_FINAL_ASSEMBLY));
2940     PetscCall(MatAssemblyEnd(is->A, MAT_FINAL_ASSEMBLY));
2941   }
2942   PetscFunctionReturn(PETSC_SUCCESS);
2943 }
2944 
2945 static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A, PetscInt n, const PetscInt rows[], PetscScalar diag, Vec x, Vec b, PetscBool columns)
2946 {
2947   Mat_IS   *matis = (Mat_IS *)A->data;
2948   PetscInt  nr, nl, len, i;
2949   PetscInt *lrows;
2950 
2951   PetscFunctionBegin;
2952   if (PetscUnlikelyDebug(columns || diag != 0. || (x && b))) {
2953     PetscBool cong;
2954 
2955     PetscCall(PetscLayoutCompare(A->rmap, A->cmap, &cong));
2956     cong = (PetscBool)(cong && matis->sf == matis->csf);
2957     PetscCheck(cong || !columns, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Columns can be zeroed if and only if A->rmap and A->cmap are congruent and the l2g maps are the same for MATIS");
2958     PetscCheck(cong || diag == 0., PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Nonzero diagonal value supported if and only if A->rmap and A->cmap are congruent and the l2g maps are the same for MATIS");
2959     PetscCheck(cong || !x || !b, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "A->rmap and A->cmap need to be congruent, and the l2g maps be the same");
2960   }
2961   /* get locally owned rows */
2962   PetscCall(PetscLayoutMapLocal(A->rmap, n, rows, &len, &lrows, NULL));
2963   /* fix right-hand side if needed */
2964   if (x && b) {
2965     const PetscScalar *xx;
2966     PetscScalar       *bb;
2967 
2968     PetscCall(VecGetArrayRead(x, &xx));
2969     PetscCall(VecGetArray(b, &bb));
2970     for (i = 0; i < len; ++i) bb[lrows[i]] = diag * xx[lrows[i]];
2971     PetscCall(VecRestoreArrayRead(x, &xx));
2972     PetscCall(VecRestoreArray(b, &bb));
2973   }
2974   /* get rows associated to the local matrices */
2975   PetscCall(MatGetSize(matis->A, &nl, NULL));
2976   PetscCall(PetscArrayzero(matis->sf_leafdata, nl));
2977   PetscCall(PetscArrayzero(matis->sf_rootdata, A->rmap->n));
2978   for (i = 0; i < len; i++) matis->sf_rootdata[lrows[i]] = 1;
2979   PetscCall(PetscFree(lrows));
2980   PetscCall(PetscSFBcastBegin(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
2981   PetscCall(PetscSFBcastEnd(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
2982   PetscCall(PetscMalloc1(nl, &lrows));
2983   for (i = 0, nr = 0; i < nl; i++)
2984     if (matis->sf_leafdata[i]) lrows[nr++] = i;
2985   PetscCall(MatISZeroRowsColumnsLocal_Private(A, nr, lrows, diag, columns));
2986   PetscCall(PetscFree(lrows));
2987   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2988   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2989   PetscFunctionReturn(PETSC_SUCCESS);
2990 }
2991 
2992 static PetscErrorCode MatZeroRows_IS(Mat A, PetscInt n, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
2993 {
2994   PetscFunctionBegin;
2995   PetscCall(MatZeroRowsColumns_Private_IS(A, n, rows, diag, x, b, PETSC_FALSE));
2996   PetscFunctionReturn(PETSC_SUCCESS);
2997 }
2998 
2999 static PetscErrorCode MatZeroRowsColumns_IS(Mat A, PetscInt n, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
3000 {
3001   PetscFunctionBegin;
3002   PetscCall(MatZeroRowsColumns_Private_IS(A, n, rows, diag, x, b, PETSC_TRUE));
3003   PetscFunctionReturn(PETSC_SUCCESS);
3004 }
3005 
3006 static PetscErrorCode MatAssemblyBegin_IS(Mat A, MatAssemblyType type)
3007 {
3008   Mat_IS *is = (Mat_IS *)A->data;
3009 
3010   PetscFunctionBegin;
3011   PetscCall(MatAssemblyBegin(is->A, type));
3012   PetscFunctionReturn(PETSC_SUCCESS);
3013 }
3014 
3015 static PetscErrorCode MatAssemblyEnd_IS(Mat A, MatAssemblyType type)
3016 {
3017   Mat_IS   *is = (Mat_IS *)A->data;
3018   PetscBool lnnz;
3019 
3020   PetscFunctionBegin;
3021   PetscCall(MatAssemblyEnd(is->A, type));
3022   /* fix for local empty rows/cols */
3023   if (is->locempty && type == MAT_FINAL_ASSEMBLY) {
3024     Mat                    newlA;
3025     ISLocalToGlobalMapping rl2g, cl2g;
3026     IS                     nzr, nzc;
3027     PetscInt               nr, nc, nnzr, nnzc;
3028     PetscBool              lnewl2g, newl2g;
3029 
3030     PetscCall(MatGetSize(is->A, &nr, &nc));
3031     PetscCall(MatFindNonzeroRowsOrCols_Basic(is->A, PETSC_FALSE, PETSC_SMALL, &nzr));
3032     if (!nzr) PetscCall(ISCreateStride(PetscObjectComm((PetscObject)is->A), nr, 0, 1, &nzr));
3033     PetscCall(MatFindNonzeroRowsOrCols_Basic(is->A, PETSC_TRUE, PETSC_SMALL, &nzc));
3034     if (!nzc) PetscCall(ISCreateStride(PetscObjectComm((PetscObject)is->A), nc, 0, 1, &nzc));
3035     PetscCall(ISGetSize(nzr, &nnzr));
3036     PetscCall(ISGetSize(nzc, &nnzc));
3037     if (nnzr != nr || nnzc != nc) { /* need new global l2g map */
3038       lnewl2g = PETSC_TRUE;
3039       PetscCallMPI(MPIU_Allreduce(&lnewl2g, &newl2g, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)A)));
3040 
3041       /* extract valid submatrix */
3042       PetscCall(MatCreateSubMatrix(is->A, nzr, nzc, MAT_INITIAL_MATRIX, &newlA));
3043     } else { /* local matrix fully populated */
3044       lnewl2g = PETSC_FALSE;
3045       PetscCallMPI(MPIU_Allreduce(&lnewl2g, &newl2g, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)A)));
3046       PetscCall(PetscObjectReference((PetscObject)is->A));
3047       newlA = is->A;
3048     }
3049 
3050     /* attach new global l2g map if needed */
3051     if (newl2g) {
3052       IS              zr, zc;
3053       const PetscInt *ridxs, *cidxs, *zridxs, *zcidxs;
3054       PetscInt       *nidxs, i;
3055 
3056       PetscCall(ISComplement(nzr, 0, nr, &zr));
3057       PetscCall(ISComplement(nzc, 0, nc, &zc));
3058       PetscCall(PetscMalloc1(PetscMax(nr, nc), &nidxs));
3059       PetscCall(ISLocalToGlobalMappingGetIndices(is->rmapping, &ridxs));
3060       PetscCall(ISLocalToGlobalMappingGetIndices(is->cmapping, &cidxs));
3061       PetscCall(ISGetIndices(zr, &zridxs));
3062       PetscCall(ISGetIndices(zc, &zcidxs));
3063       PetscCall(ISGetLocalSize(zr, &nnzr));
3064       PetscCall(ISGetLocalSize(zc, &nnzc));
3065 
3066       PetscCall(PetscArraycpy(nidxs, ridxs, nr));
3067       for (i = 0; i < nnzr; i++) nidxs[zridxs[i]] = -1;
3068       PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A), 1, nr, nidxs, PETSC_COPY_VALUES, &rl2g));
3069       PetscCall(PetscArraycpy(nidxs, cidxs, nc));
3070       for (i = 0; i < nnzc; i++) nidxs[zcidxs[i]] = -1;
3071       PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A), 1, nc, nidxs, PETSC_COPY_VALUES, &cl2g));
3072 
3073       PetscCall(ISRestoreIndices(zr, &zridxs));
3074       PetscCall(ISRestoreIndices(zc, &zcidxs));
3075       PetscCall(ISLocalToGlobalMappingRestoreIndices(is->rmapping, &ridxs));
3076       PetscCall(ISLocalToGlobalMappingRestoreIndices(is->cmapping, &cidxs));
3077       PetscCall(ISDestroy(&nzr));
3078       PetscCall(ISDestroy(&nzc));
3079       PetscCall(ISDestroy(&zr));
3080       PetscCall(ISDestroy(&zc));
3081       PetscCall(PetscFree(nidxs));
3082       PetscCall(MatSetLocalToGlobalMapping(A, rl2g, cl2g));
3083       PetscCall(ISLocalToGlobalMappingDestroy(&rl2g));
3084       PetscCall(ISLocalToGlobalMappingDestroy(&cl2g));
3085     }
3086     PetscCall(MatISSetLocalMat(A, newlA));
3087     PetscCall(MatDestroy(&newlA));
3088     PetscCall(ISDestroy(&nzr));
3089     PetscCall(ISDestroy(&nzc));
3090     is->locempty = PETSC_FALSE;
3091   }
3092   lnnz          = (PetscBool)(is->A->nonzerostate == is->lnnzstate);
3093   is->lnnzstate = is->A->nonzerostate;
3094   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &lnnz, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
3095   if (!lnnz) A->nonzerostate++;
3096   PetscFunctionReturn(PETSC_SUCCESS);
3097 }
3098 
3099 static PetscErrorCode MatISGetLocalMat_IS(Mat mat, Mat *local)
3100 {
3101   Mat_IS *is = (Mat_IS *)mat->data;
3102 
3103   PetscFunctionBegin;
3104   *local = is->A;
3105   PetscFunctionReturn(PETSC_SUCCESS);
3106 }
3107 
3108 static PetscErrorCode MatISRestoreLocalMat_IS(Mat mat, Mat *local)
3109 {
3110   PetscFunctionBegin;
3111   *local = NULL;
3112   PetscFunctionReturn(PETSC_SUCCESS);
3113 }
3114 
3115 /*@
3116   MatISGetLocalMat - Gets the local matrix stored inside a `MATIS` matrix.
3117 
3118   Not Collective.
3119 
3120   Input Parameter:
3121 . mat - the matrix
3122 
3123   Output Parameter:
3124 . local - the local matrix
3125 
3126   Level: intermediate
3127 
3128   Notes:
3129   This can be called if you have precomputed the nonzero structure of the
3130   matrix and want to provide it to the inner matrix object to improve the performance
3131   of the `MatSetValues()` operation.
3132 
3133   Call `MatISRestoreLocalMat()` when finished with the local matrix.
3134 
3135 .seealso: [](ch_matrices), `Mat`, `MATIS`, `MatISRestoreLocalMat()`
3136 @*/
3137 PetscErrorCode MatISGetLocalMat(Mat mat, Mat *local)
3138 {
3139   PetscFunctionBegin;
3140   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
3141   PetscAssertPointer(local, 2);
3142   PetscUseMethod(mat, "MatISGetLocalMat_C", (Mat, Mat *), (mat, local));
3143   PetscFunctionReturn(PETSC_SUCCESS);
3144 }
3145 
3146 /*@
3147   MatISRestoreLocalMat - Restores the local matrix obtained with `MatISGetLocalMat()`
3148 
3149   Not Collective.
3150 
3151   Input Parameters:
3152 + mat   - the matrix
3153 - local - the local matrix
3154 
3155   Level: intermediate
3156 
3157 .seealso: [](ch_matrices), `Mat`, `MATIS`, `MatISGetLocalMat()`
3158 @*/
3159 PetscErrorCode MatISRestoreLocalMat(Mat mat, Mat *local)
3160 {
3161   PetscFunctionBegin;
3162   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
3163   PetscAssertPointer(local, 2);
3164   PetscUseMethod(mat, "MatISRestoreLocalMat_C", (Mat, Mat *), (mat, local));
3165   PetscFunctionReturn(PETSC_SUCCESS);
3166 }
3167 
3168 static PetscErrorCode MatISSetLocalMatType_IS(Mat mat, MatType mtype)
3169 {
3170   Mat_IS *is = (Mat_IS *)mat->data;
3171 
3172   PetscFunctionBegin;
3173   if (is->A) PetscCall(MatSetType(is->A, mtype));
3174   PetscCall(PetscFree(is->lmattype));
3175   PetscCall(PetscStrallocpy(mtype, &is->lmattype));
3176   PetscFunctionReturn(PETSC_SUCCESS);
3177 }
3178 
3179 /*@
3180   MatISSetLocalMatType - Specifies the type of local matrix inside the `MATIS`
3181 
3182   Logically Collective.
3183 
3184   Input Parameters:
3185 + mat   - the matrix
3186 - mtype - the local matrix type
3187 
3188   Level: intermediate
3189 
3190 .seealso: [](ch_matrices), `Mat`, `MATIS`, `MatSetType()`, `MatType`
3191 @*/
3192 PetscErrorCode MatISSetLocalMatType(Mat mat, MatType mtype)
3193 {
3194   PetscFunctionBegin;
3195   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
3196   PetscUseMethod(mat, "MatISSetLocalMatType_C", (Mat, MatType), (mat, mtype));
3197   PetscFunctionReturn(PETSC_SUCCESS);
3198 }
3199 
3200 static PetscErrorCode MatISSetLocalMat_IS(Mat mat, Mat local)
3201 {
3202   Mat_IS   *is = (Mat_IS *)mat->data;
3203   PetscInt  nrows, ncols, orows, ocols;
3204   MatType   mtype, otype;
3205   PetscBool sametype = PETSC_TRUE;
3206 
3207   PetscFunctionBegin;
3208   if (is->A && !is->islocalref) {
3209     PetscCall(MatGetSize(is->A, &orows, &ocols));
3210     PetscCall(MatGetSize(local, &nrows, &ncols));
3211     PetscCheck(orows == nrows && ocols == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local MATIS matrix should be of size %" PetscInt_FMT "x%" PetscInt_FMT " (passed a %" PetscInt_FMT "x%" PetscInt_FMT " matrix)", orows, ocols, nrows, ncols);
3212     PetscCall(MatGetType(local, &mtype));
3213     PetscCall(MatGetType(is->A, &otype));
3214     PetscCall(PetscStrcmp(mtype, otype, &sametype));
3215   }
3216   PetscCall(PetscObjectReference((PetscObject)local));
3217   PetscCall(MatDestroy(&is->A));
3218   is->A = local;
3219   PetscCall(MatGetType(is->A, &mtype));
3220   PetscCall(MatISSetLocalMatType(mat, mtype));
3221   if (!sametype && !is->islocalref) PetscCall(MatISSetUpScatters_Private(mat));
3222   is->lnnzstate = 0;
3223   PetscFunctionReturn(PETSC_SUCCESS);
3224 }
3225 
3226 /*@
3227   MatISSetLocalMat - Replace the local matrix stored inside a `MATIS` object.
3228 
3229   Not Collective
3230 
3231   Input Parameters:
3232 + mat   - the matrix
3233 - local - the local matrix
3234 
3235   Level: intermediate
3236 
3237 .seealso: [](ch_matrices), `Mat`, `MATIS`, `MatISSetLocalMatType`, `MatISGetLocalMat()`
3238 @*/
3239 PetscErrorCode MatISSetLocalMat(Mat mat, Mat local)
3240 {
3241   PetscFunctionBegin;
3242   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
3243   PetscValidHeaderSpecific(local, MAT_CLASSID, 2);
3244   PetscUseMethod(mat, "MatISSetLocalMat_C", (Mat, Mat), (mat, local));
3245   PetscFunctionReturn(PETSC_SUCCESS);
3246 }
3247 
3248 static PetscErrorCode MatZeroEntries_IS(Mat A)
3249 {
3250   Mat_IS *a = (Mat_IS *)A->data;
3251 
3252   PetscFunctionBegin;
3253   PetscCall(MatZeroEntries(a->A));
3254   PetscFunctionReturn(PETSC_SUCCESS);
3255 }
3256 
3257 static PetscErrorCode MatScale_IS(Mat A, PetscScalar a)
3258 {
3259   Mat_IS *is = (Mat_IS *)A->data;
3260 
3261   PetscFunctionBegin;
3262   PetscCall(MatScale(is->A, a));
3263   PetscFunctionReturn(PETSC_SUCCESS);
3264 }
3265 
3266 static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
3267 {
3268   Mat_IS *is = (Mat_IS *)A->data;
3269 
3270   PetscFunctionBegin;
3271   /* get diagonal of the local matrix */
3272   PetscCall(MatGetDiagonal(is->A, is->y));
3273 
3274   /* scatter diagonal back into global vector */
3275   PetscCall(VecSet(v, 0));
3276   PetscCall(VecScatterBegin(is->rctx, is->y, v, ADD_VALUES, SCATTER_REVERSE));
3277   PetscCall(VecScatterEnd(is->rctx, is->y, v, ADD_VALUES, SCATTER_REVERSE));
3278   PetscFunctionReturn(PETSC_SUCCESS);
3279 }
3280 
3281 static PetscErrorCode MatSetOption_IS(Mat A, MatOption op, PetscBool flg)
3282 {
3283   Mat_IS *a = (Mat_IS *)A->data;
3284 
3285   PetscFunctionBegin;
3286   PetscCall(MatSetOption(a->A, op, flg));
3287   PetscFunctionReturn(PETSC_SUCCESS);
3288 }
3289 
3290 static PetscErrorCode MatAXPY_IS(Mat Y, PetscScalar a, Mat X, MatStructure str)
3291 {
3292   Mat_IS *y = (Mat_IS *)Y->data;
3293   Mat_IS *x;
3294 
3295   PetscFunctionBegin;
3296   if (PetscDefined(USE_DEBUG)) {
3297     PetscBool ismatis;
3298     PetscCall(PetscObjectTypeCompare((PetscObject)X, MATIS, &ismatis));
3299     PetscCheck(ismatis, PetscObjectComm((PetscObject)Y), PETSC_ERR_SUP, "Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS");
3300   }
3301   x = (Mat_IS *)X->data;
3302   PetscCall(MatAXPY(y->A, a, x->A, str));
3303   PetscFunctionReturn(PETSC_SUCCESS);
3304 }
3305 
3306 static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A, IS row, IS col, Mat *submat)
3307 {
3308   Mat                    lA;
3309   Mat_IS                *matis = (Mat_IS *)A->data;
3310   ISLocalToGlobalMapping rl2g, cl2g;
3311   IS                     is;
3312   const PetscInt        *rg, *rl;
3313   PetscInt               nrg;
3314   PetscInt               N, M, nrl, i, *idxs;
3315 
3316   PetscFunctionBegin;
3317   PetscCall(ISLocalToGlobalMappingGetIndices(A->rmap->mapping, &rg));
3318   PetscCall(ISGetLocalSize(row, &nrl));
3319   PetscCall(ISGetIndices(row, &rl));
3320   PetscCall(ISLocalToGlobalMappingGetSize(A->rmap->mapping, &nrg));
3321   if (PetscDefined(USE_DEBUG)) {
3322     for (i = 0; i < nrl; i++) PetscCheck(rl[i] < nrg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Local row index %" PetscInt_FMT " -> %" PetscInt_FMT " greater then maximum possible %" PetscInt_FMT, i, rl[i], nrg);
3323   }
3324   PetscCall(PetscMalloc1(nrg, &idxs));
3325   /* map from [0,nrl) to row */
3326   for (i = 0; i < nrl; i++) idxs[i] = rl[i];
3327   for (i = nrl; i < nrg; i++) idxs[i] = -1;
3328   PetscCall(ISRestoreIndices(row, &rl));
3329   PetscCall(ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping, &rg));
3330   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), nrg, idxs, PETSC_OWN_POINTER, &is));
3331   PetscCall(ISLocalToGlobalMappingCreateIS(is, &rl2g));
3332   PetscCall(ISDestroy(&is));
3333   /* compute new l2g map for columns */
3334   if (col != row || matis->rmapping != matis->cmapping || matis->A->rmap->mapping != matis->A->cmap->mapping) {
3335     const PetscInt *cg, *cl;
3336     PetscInt        ncg;
3337     PetscInt        ncl;
3338 
3339     PetscCall(ISLocalToGlobalMappingGetIndices(A->cmap->mapping, &cg));
3340     PetscCall(ISGetLocalSize(col, &ncl));
3341     PetscCall(ISGetIndices(col, &cl));
3342     PetscCall(ISLocalToGlobalMappingGetSize(A->cmap->mapping, &ncg));
3343     if (PetscDefined(USE_DEBUG)) {
3344       for (i = 0; i < ncl; i++) PetscCheck(cl[i] < ncg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Local column index %" PetscInt_FMT " -> %" PetscInt_FMT " greater then maximum possible %" PetscInt_FMT, i, cl[i], ncg);
3345     }
3346     PetscCall(PetscMalloc1(ncg, &idxs));
3347     /* map from [0,ncl) to col */
3348     for (i = 0; i < ncl; i++) idxs[i] = cl[i];
3349     for (i = ncl; i < ncg; i++) idxs[i] = -1;
3350     PetscCall(ISRestoreIndices(col, &cl));
3351     PetscCall(ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping, &cg));
3352     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), ncg, idxs, PETSC_OWN_POINTER, &is));
3353     PetscCall(ISLocalToGlobalMappingCreateIS(is, &cl2g));
3354     PetscCall(ISDestroy(&is));
3355   } else {
3356     PetscCall(PetscObjectReference((PetscObject)rl2g));
3357     cl2g = rl2g;
3358   }
3359   /* create the MATIS submatrix */
3360   PetscCall(MatGetSize(A, &M, &N));
3361   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), submat));
3362   PetscCall(MatSetSizes(*submat, PETSC_DECIDE, PETSC_DECIDE, M, N));
3363   PetscCall(MatSetType(*submat, MATIS));
3364   matis             = (Mat_IS *)((*submat)->data);
3365   matis->islocalref = PETSC_TRUE;
3366   PetscCall(MatSetLocalToGlobalMapping(*submat, rl2g, cl2g));
3367   PetscCall(MatISGetLocalMat(A, &lA));
3368   PetscCall(MatISSetLocalMat(*submat, lA));
3369   PetscCall(MatSetUp(*submat));
3370   PetscCall(MatAssemblyBegin(*submat, MAT_FINAL_ASSEMBLY));
3371   PetscCall(MatAssemblyEnd(*submat, MAT_FINAL_ASSEMBLY));
3372   PetscCall(ISLocalToGlobalMappingDestroy(&rl2g));
3373   PetscCall(ISLocalToGlobalMappingDestroy(&cl2g));
3374 
3375   /* remove unsupported ops */
3376   PetscCall(PetscMemzero((*submat)->ops, sizeof(struct _MatOps)));
3377   (*submat)->ops->destroy               = MatDestroy_IS;
3378   (*submat)->ops->setvalueslocal        = MatSetValuesLocal_SubMat_IS;
3379   (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS;
3380   (*submat)->ops->assemblybegin         = MatAssemblyBegin_IS;
3381   (*submat)->ops->assemblyend           = MatAssemblyEnd_IS;
3382   PetscFunctionReturn(PETSC_SUCCESS);
3383 }
3384 
3385 static PetscErrorCode MatSetFromOptions_IS(Mat A, PetscOptionItems PetscOptionsObject)
3386 {
3387   Mat_IS   *a = (Mat_IS *)A->data;
3388   char      type[256];
3389   PetscBool flg;
3390 
3391   PetscFunctionBegin;
3392   PetscOptionsHeadBegin(PetscOptionsObject, "MATIS options");
3393   PetscCall(PetscOptionsDeprecated("-matis_keepassembled", "-mat_is_keepassembled", "3.21", NULL));
3394   PetscCall(PetscOptionsDeprecated("-matis_fixempty", "-mat_is_fixempty", "3.21", NULL));
3395   PetscCall(PetscOptionsDeprecated("-matis_storel2l", "-mat_is_storel2l", "3.21", NULL));
3396   PetscCall(PetscOptionsDeprecated("-matis_localmat_type", "-mat_is_localmat_type", "3.21", NULL));
3397   PetscCall(PetscOptionsBool("-mat_is_keepassembled", "Store an assembled version if needed", NULL, a->keepassembled, &a->keepassembled, NULL));
3398   PetscCall(PetscOptionsBool("-mat_is_fixempty", "Fix local matrices in case of empty local rows/columns", "MatISFixLocalEmpty", a->locempty, &a->locempty, NULL));
3399   PetscCall(PetscOptionsBool("-mat_is_storel2l", "Store local-to-local matrices generated from PtAP operations", "MatISStoreL2L", a->storel2l, &a->storel2l, NULL));
3400   PetscCall(PetscOptionsBool("-mat_is_allow_repeated", "Allow local repeated entries", "MatISSetAllowRepeated", a->allow_repeated, &a->allow_repeated, NULL));
3401   PetscCall(PetscOptionsFList("-mat_is_localmat_type", "Matrix type", "MatISSetLocalMatType", MatList, a->lmattype, type, 256, &flg));
3402   if (flg) PetscCall(MatISSetLocalMatType(A, type));
3403   if (a->A) PetscCall(MatSetFromOptions(a->A));
3404   PetscOptionsHeadEnd();
3405   PetscFunctionReturn(PETSC_SUCCESS);
3406 }
3407 
3408 /*@
3409   MatCreateIS - Creates a "process" unassembled matrix.
3410 
3411   Collective.
3412 
3413   Input Parameters:
3414 + comm - MPI communicator that will share the matrix
3415 . bs   - block size of the matrix
3416 . m    - local size of left vector used in matrix vector products
3417 . n    - local size of right vector used in matrix vector products
3418 . M    - global size of left vector used in matrix vector products
3419 . N    - global size of right vector used in matrix vector products
3420 . rmap - local to global map for rows
3421 - cmap - local to global map for cols
3422 
3423   Output Parameter:
3424 . A - the resulting matrix
3425 
3426   Level: intermediate
3427 
3428   Notes:
3429   `m` and `n` are NOT related to the size of the map; they represent the size of the local parts of the distributed vectors
3430   used in `MatMult()` operations. The local sizes of `rmap` and `cmap` define the size of the local matrices.
3431 
3432   If `rmap` (`cmap`) is `NULL`, then the local row (column) spaces matches the global space.
3433 
3434 .seealso: [](ch_matrices), `Mat`, `MATIS`, `MatSetLocalToGlobalMapping()`
3435 @*/
3436 PetscErrorCode MatCreateIS(MPI_Comm comm, PetscInt bs, PetscInt m, PetscInt n, PetscInt M, PetscInt N, ISLocalToGlobalMapping rmap, ISLocalToGlobalMapping cmap, Mat *A)
3437 {
3438   PetscFunctionBegin;
3439   PetscCall(MatCreate(comm, A));
3440   PetscCall(MatSetSizes(*A, m, n, M, N));
3441   if (bs > 0) PetscCall(MatSetBlockSize(*A, bs));
3442   PetscCall(MatSetType(*A, MATIS));
3443   PetscCall(MatSetLocalToGlobalMapping(*A, rmap, cmap));
3444   PetscFunctionReturn(PETSC_SUCCESS);
3445 }
3446 
3447 static PetscErrorCode MatHasOperation_IS(Mat A, MatOperation op, PetscBool *has)
3448 {
3449   Mat_IS      *a              = (Mat_IS *)A->data;
3450   MatOperation tobefiltered[] = {MATOP_MULT_ADD, MATOP_MULT_TRANSPOSE_ADD, MATOP_GET_DIAGONAL_BLOCK, MATOP_INCREASE_OVERLAP};
3451 
3452   PetscFunctionBegin;
3453   *has = PETSC_FALSE;
3454   if (!((void **)A->ops)[op] || !a->A) PetscFunctionReturn(PETSC_SUCCESS);
3455   *has = PETSC_TRUE;
3456   for (PetscInt i = 0; i < (PetscInt)PETSC_STATIC_ARRAY_LENGTH(tobefiltered); i++)
3457     if (op == tobefiltered[i]) PetscFunctionReturn(PETSC_SUCCESS);
3458   PetscCall(MatHasOperation(a->A, op, has));
3459   PetscFunctionReturn(PETSC_SUCCESS);
3460 }
3461 
3462 static PetscErrorCode MatSetValuesCOO_IS(Mat A, const PetscScalar v[], InsertMode imode)
3463 {
3464   Mat_IS *a = (Mat_IS *)A->data;
3465 
3466   PetscFunctionBegin;
3467   PetscCall(MatSetValuesCOO(a->A, v, imode));
3468   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3469   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
3470   PetscFunctionReturn(PETSC_SUCCESS);
3471 }
3472 
3473 static PetscErrorCode MatSetPreallocationCOOLocal_IS(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
3474 {
3475   Mat_IS *a = (Mat_IS *)A->data;
3476 
3477   PetscFunctionBegin;
3478   PetscCheck(a->A, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to provide l2g map first via MatSetLocalToGlobalMapping");
3479   if (a->A->rmap->mapping || a->A->cmap->mapping) {
3480     PetscCall(MatSetPreallocationCOOLocal(a->A, ncoo, coo_i, coo_j));
3481   } else {
3482     PetscCall(MatSetPreallocationCOO(a->A, ncoo, coo_i, coo_j));
3483   }
3484   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetValuesCOO_C", MatSetValuesCOO_IS));
3485   A->preallocated = PETSC_TRUE;
3486   PetscFunctionReturn(PETSC_SUCCESS);
3487 }
3488 
3489 static PetscErrorCode MatSetPreallocationCOO_IS(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
3490 {
3491   Mat_IS  *a = (Mat_IS *)A->data;
3492   PetscInt ncoo_i;
3493 
3494   PetscFunctionBegin;
3495   PetscCheck(a->A, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Need to provide l2g map first via MatSetLocalToGlobalMapping");
3496   PetscCall(PetscIntCast(ncoo, &ncoo_i));
3497   PetscCall(ISGlobalToLocalMappingApply(a->rmapping, IS_GTOLM_MASK, ncoo_i, coo_i, NULL, coo_i));
3498   PetscCall(ISGlobalToLocalMappingApply(a->cmapping, IS_GTOLM_MASK, ncoo_i, coo_j, NULL, coo_j));
3499   PetscCall(MatSetPreallocationCOO(a->A, ncoo, coo_i, coo_j));
3500   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetValuesCOO_C", MatSetValuesCOO_IS));
3501   A->preallocated = PETSC_TRUE;
3502   PetscFunctionReturn(PETSC_SUCCESS);
3503 }
3504 
3505 static PetscErrorCode MatISGetAssembled_Private(Mat A, Mat *tA)
3506 {
3507   Mat_IS          *a = (Mat_IS *)A->data;
3508   PetscObjectState Astate, aAstate       = PETSC_INT_MIN;
3509   PetscObjectState Annzstate, aAnnzstate = PETSC_INT_MIN;
3510 
3511   PetscFunctionBegin;
3512   PetscCall(PetscObjectStateGet((PetscObject)A, &Astate));
3513   Annzstate = A->nonzerostate;
3514   if (a->assembledA) {
3515     PetscCall(PetscObjectStateGet((PetscObject)a->assembledA, &aAstate));
3516     aAnnzstate = a->assembledA->nonzerostate;
3517   }
3518   if (aAnnzstate != Annzstate) PetscCall(MatDestroy(&a->assembledA));
3519   if (Astate != aAstate || !a->assembledA) {
3520     MatType     aAtype;
3521     PetscMPIInt size;
3522     PetscInt    rbs, cbs, bs;
3523 
3524     /* the assembled form is used as temporary storage for parallel operations
3525        like createsubmatrices and the like, do not waste device memory */
3526     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
3527     PetscCall(ISLocalToGlobalMappingGetBlockSize(a->cmapping, &cbs));
3528     PetscCall(ISLocalToGlobalMappingGetBlockSize(a->rmapping, &rbs));
3529     bs = rbs == cbs ? rbs : 1;
3530     if (a->assembledA) PetscCall(MatGetType(a->assembledA, &aAtype));
3531     else if (size > 1) aAtype = bs > 1 ? MATMPIBAIJ : MATMPIAIJ;
3532     else aAtype = bs > 1 ? MATSEQBAIJ : MATSEQAIJ;
3533 
3534     PetscCall(MatConvert(A, aAtype, a->assembledA ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX, &a->assembledA));
3535     PetscCall(PetscObjectStateSet((PetscObject)a->assembledA, Astate));
3536     a->assembledA->nonzerostate = Annzstate;
3537   }
3538   PetscCall(PetscObjectReference((PetscObject)a->assembledA));
3539   *tA = a->assembledA;
3540   if (!a->keepassembled) PetscCall(MatDestroy(&a->assembledA));
3541   PetscFunctionReturn(PETSC_SUCCESS);
3542 }
3543 
3544 static PetscErrorCode MatISRestoreAssembled_Private(Mat A, Mat *tA)
3545 {
3546   PetscFunctionBegin;
3547   PetscCall(MatDestroy(tA));
3548   *tA = NULL;
3549   PetscFunctionReturn(PETSC_SUCCESS);
3550 }
3551 
3552 static PetscErrorCode MatGetDiagonalBlock_IS(Mat A, Mat *dA)
3553 {
3554   Mat_IS          *a = (Mat_IS *)A->data;
3555   PetscObjectState Astate, dAstate = PETSC_INT_MIN;
3556 
3557   PetscFunctionBegin;
3558   PetscCall(PetscObjectStateGet((PetscObject)A, &Astate));
3559   if (a->dA) PetscCall(PetscObjectStateGet((PetscObject)a->dA, &dAstate));
3560   if (Astate != dAstate) {
3561     Mat     tA;
3562     MatType ltype;
3563 
3564     PetscCall(MatDestroy(&a->dA));
3565     PetscCall(MatISGetAssembled_Private(A, &tA));
3566     PetscCall(MatGetDiagonalBlock(tA, &a->dA));
3567     PetscCall(MatPropagateSymmetryOptions(tA, a->dA));
3568     PetscCall(MatGetType(a->A, &ltype));
3569     PetscCall(MatConvert(a->dA, ltype, MAT_INPLACE_MATRIX, &a->dA));
3570     PetscCall(PetscObjectReference((PetscObject)a->dA));
3571     PetscCall(MatISRestoreAssembled_Private(A, &tA));
3572     PetscCall(PetscObjectStateSet((PetscObject)a->dA, Astate));
3573   }
3574   *dA = a->dA;
3575   PetscFunctionReturn(PETSC_SUCCESS);
3576 }
3577 
3578 static PetscErrorCode MatCreateSubMatrices_IS(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse reuse, Mat *submat[])
3579 {
3580   Mat tA;
3581 
3582   PetscFunctionBegin;
3583   PetscCall(MatISGetAssembled_Private(A, &tA));
3584   PetscCall(MatCreateSubMatrices(tA, n, irow, icol, reuse, submat));
3585   /* MatCreateSubMatrices_MPIAIJ is a mess at the moment */
3586 #if 0
3587   {
3588     Mat_IS    *a = (Mat_IS*)A->data;
3589     MatType   ltype;
3590     VecType   vtype;
3591     char      *flg;
3592 
3593     PetscCall(MatGetType(a->A,&ltype));
3594     PetscCall(MatGetVecType(a->A,&vtype));
3595     PetscCall(PetscStrstr(vtype,"cuda",&flg));
3596     if (!flg) PetscCall(PetscStrstr(vtype,"hip",&flg));
3597     if (!flg) PetscCall(PetscStrstr(vtype,"kokkos",&flg));
3598     if (flg) {
3599       for (PetscInt i = 0; i < n; i++) {
3600         Mat sA = (*submat)[i];
3601 
3602         PetscCall(MatConvert(sA,ltype,MAT_INPLACE_MATRIX,&sA));
3603         (*submat)[i] = sA;
3604       }
3605     }
3606   }
3607 #endif
3608   PetscCall(MatISRestoreAssembled_Private(A, &tA));
3609   PetscFunctionReturn(PETSC_SUCCESS);
3610 }
3611 
3612 static PetscErrorCode MatIncreaseOverlap_IS(Mat A, PetscInt n, IS is[], PetscInt ov)
3613 {
3614   Mat tA;
3615 
3616   PetscFunctionBegin;
3617   PetscCall(MatISGetAssembled_Private(A, &tA));
3618   PetscCall(MatIncreaseOverlap(tA, n, is, ov));
3619   PetscCall(MatISRestoreAssembled_Private(A, &tA));
3620   PetscFunctionReturn(PETSC_SUCCESS);
3621 }
3622 
3623 /*@
3624   MatISGetLocalToGlobalMapping - Gets the local-to-global numbering of the `MATIS` object
3625 
3626   Not Collective
3627 
3628   Input Parameter:
3629 . A - the matrix
3630 
3631   Output Parameters:
3632 + rmapping - row mapping
3633 - cmapping - column mapping
3634 
3635   Level: advanced
3636 
3637   Note:
3638   The returned map can be different from the one used to construct the `MATIS` object, since it will not contain negative or repeated indices.
3639 
3640 .seealso: [](ch_matrices), `Mat`, `MATIS`, `MatSetLocalToGlobalMapping()`
3641 @*/
3642 PetscErrorCode MatISGetLocalToGlobalMapping(Mat A, ISLocalToGlobalMapping *rmapping, ISLocalToGlobalMapping *cmapping)
3643 {
3644   PetscFunctionBegin;
3645   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
3646   PetscValidType(A, 1);
3647   if (rmapping) PetscAssertPointer(rmapping, 2);
3648   if (cmapping) PetscAssertPointer(cmapping, 3);
3649   PetscUseMethod(A, "MatISGetLocalToGlobalMapping_C", (Mat, ISLocalToGlobalMapping *, ISLocalToGlobalMapping *), (A, rmapping, cmapping));
3650   PetscFunctionReturn(PETSC_SUCCESS);
3651 }
3652 
3653 static PetscErrorCode MatISGetLocalToGlobalMapping_IS(Mat A, ISLocalToGlobalMapping *r, ISLocalToGlobalMapping *c)
3654 {
3655   Mat_IS *a = (Mat_IS *)A->data;
3656 
3657   PetscFunctionBegin;
3658   if (r) *r = a->rmapping;
3659   if (c) *c = a->cmapping;
3660   PetscFunctionReturn(PETSC_SUCCESS);
3661 }
3662 
3663 static PetscErrorCode MatSetBlockSizes_IS(Mat A, PetscInt rbs, PetscInt cbs)
3664 {
3665   Mat_IS *a = (Mat_IS *)A->data;
3666 
3667   PetscFunctionBegin;
3668   if (a->A) PetscCall(MatSetBlockSizes(a->A, rbs, cbs));
3669   PetscFunctionReturn(PETSC_SUCCESS);
3670 }
3671 
3672 /*MC
3673   MATIS - MATIS = "is" - A matrix type to be used for non-overlapping domain decomposition methods (e.g. `PCBDDC` or `KSPFETIDP`).
3674   This stores the matrices in globally unassembled form and the parallel matrix vector product is handled "implicitly".
3675 
3676   Options Database Keys:
3677 + -mat_type is           - Set the matrix type to `MATIS`.
3678 . -mat_is_allow_repeated - Allow repeated entries in the local part of the local to global maps.
3679 . -mat_is_fixempty       - Fix local matrices in case of empty local rows/columns.
3680 - -mat_is_storel2l       - Store the local-to-local operators generated by the Galerkin process of `MatPtAP()`.
3681 
3682   Level: intermediate
3683 
3684   Notes:
3685   Options prefix for the inner matrix are given by `-is_mat_xxx`
3686 
3687   You must call `MatSetLocalToGlobalMapping()` before using this matrix type.
3688 
3689   You can do matrix preallocation on the local matrix after you obtain it with
3690   `MatISGetLocalMat()`; otherwise, you could use `MatISSetPreallocation()` or `MatXAIJSetPreallocation()`
3691 
3692 .seealso: [](ch_matrices), `Mat`, `MATIS`, `Mat`, `MatISGetLocalMat()`, `MatSetLocalToGlobalMapping()`, `MatISSetPreallocation()`, `MatCreateIS()`, `PCBDDC`, `KSPFETIDP`
3693 M*/
3694 PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
3695 {
3696   Mat_IS *a;
3697 
3698   PetscFunctionBegin;
3699   PetscCall(PetscNew(&a));
3700   PetscCall(PetscStrallocpy(MATAIJ, &a->lmattype));
3701   A->data = (void *)a;
3702 
3703   /* matrix ops */
3704   PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps)));
3705   A->ops->mult                    = MatMult_IS;
3706   A->ops->multadd                 = MatMultAdd_IS;
3707   A->ops->multtranspose           = MatMultTranspose_IS;
3708   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
3709   A->ops->destroy                 = MatDestroy_IS;
3710   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
3711   A->ops->setvalues               = MatSetValues_IS;
3712   A->ops->setvaluesblocked        = MatSetValuesBlocked_IS;
3713   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
3714   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
3715   A->ops->zerorows                = MatZeroRows_IS;
3716   A->ops->zerorowscolumns         = MatZeroRowsColumns_IS;
3717   A->ops->assemblybegin           = MatAssemblyBegin_IS;
3718   A->ops->assemblyend             = MatAssemblyEnd_IS;
3719   A->ops->view                    = MatView_IS;
3720   A->ops->load                    = MatLoad_IS;
3721   A->ops->zeroentries             = MatZeroEntries_IS;
3722   A->ops->scale                   = MatScale_IS;
3723   A->ops->getdiagonal             = MatGetDiagonal_IS;
3724   A->ops->setoption               = MatSetOption_IS;
3725   A->ops->ishermitian             = MatIsHermitian_IS;
3726   A->ops->issymmetric             = MatIsSymmetric_IS;
3727   A->ops->isstructurallysymmetric = MatIsStructurallySymmetric_IS;
3728   A->ops->duplicate               = MatDuplicate_IS;
3729   A->ops->missingdiagonal         = MatMissingDiagonal_IS;
3730   A->ops->copy                    = MatCopy_IS;
3731   A->ops->getlocalsubmatrix       = MatGetLocalSubMatrix_IS;
3732   A->ops->createsubmatrix         = MatCreateSubMatrix_IS;
3733   A->ops->axpy                    = MatAXPY_IS;
3734   A->ops->diagonalset             = MatDiagonalSet_IS;
3735   A->ops->shift                   = MatShift_IS;
3736   A->ops->transpose               = MatTranspose_IS;
3737   A->ops->getinfo                 = MatGetInfo_IS;
3738   A->ops->diagonalscale           = MatDiagonalScale_IS;
3739   A->ops->setfromoptions          = MatSetFromOptions_IS;
3740   A->ops->setup                   = MatSetUp_IS;
3741   A->ops->hasoperation            = MatHasOperation_IS;
3742   A->ops->getdiagonalblock        = MatGetDiagonalBlock_IS;
3743   A->ops->createsubmatrices       = MatCreateSubMatrices_IS;
3744   A->ops->increaseoverlap         = MatIncreaseOverlap_IS;
3745   A->ops->setblocksizes           = MatSetBlockSizes_IS;
3746 
3747   /* special MATIS functions */
3748   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISSetLocalMatType_C", MatISSetLocalMatType_IS));
3749   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISGetLocalMat_C", MatISGetLocalMat_IS));
3750   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISRestoreLocalMat_C", MatISRestoreLocalMat_IS));
3751   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISSetLocalMat_C", MatISSetLocalMat_IS));
3752   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISSetPreallocation_C", MatISSetPreallocation_IS));
3753   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISSetAllowRepeated_C", MatISSetAllowRepeated_IS));
3754   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISStoreL2L_C", MatISStoreL2L_IS));
3755   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISFixLocalEmpty_C", MatISFixLocalEmpty_IS));
3756   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatISGetLocalToGlobalMapping_C", MatISGetLocalToGlobalMapping_IS));
3757   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_mpiaij_C", MatConvert_IS_XAIJ));
3758   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_mpibaij_C", MatConvert_IS_XAIJ));
3759   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_mpisbaij_C", MatConvert_IS_XAIJ));
3760   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_seqaij_C", MatConvert_IS_XAIJ));
3761   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_seqbaij_C", MatConvert_IS_XAIJ));
3762   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_seqsbaij_C", MatConvert_IS_XAIJ));
3763   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_is_aij_C", MatConvert_IS_XAIJ));
3764   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", MatSetPreallocationCOOLocal_IS));
3765   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetPreallocationCOO_C", MatSetPreallocationCOO_IS));
3766   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATIS));
3767   PetscFunctionReturn(PETSC_SUCCESS);
3768 }
3769