xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 1fcf445a439848b2e538efa39b3f63deb630163e)
1 /*
2     Defines the basic matrix operations for the AIJ (compressed row)
3   matrix storage format.
4 */
5 
6 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
7 #include <petscblaslapack.h>
8 #include <petscbt.h>
9 #include <petsc/private/kernels/blocktranspose.h>
10 
11 /* defines MatSetValues_Seq_Hash(), MatAssemblyEnd_Seq_Hash(), MatSetUp_Seq_Hash() */
12 #define TYPE AIJ
13 #define TYPE_BS
14 #include "../src/mat/impls/aij/seq/seqhashmatsetvalues.h"
15 #include "../src/mat/impls/aij/seq/seqhashmat.h"
16 #undef TYPE
17 #undef TYPE_BS
18 
19 PetscErrorCode MatSeqAIJSetTypeFromOptions(Mat A)
20 {
21   PetscBool flg;
22   char      type[256];
23 
24   PetscFunctionBegin;
25   PetscObjectOptionsBegin((PetscObject)A);
26   PetscCall(PetscOptionsFList("-mat_seqaij_type", "Matrix SeqAIJ type", "MatSeqAIJSetType", MatSeqAIJList, "seqaij", type, 256, &flg));
27   if (flg) PetscCall(MatSeqAIJSetType(A, type));
28   PetscOptionsEnd();
29   PetscFunctionReturn(PETSC_SUCCESS);
30 }
31 
32 PetscErrorCode MatGetColumnReductions_SeqAIJ(Mat A, PetscInt type, PetscReal *reductions)
33 {
34   PetscInt    i, m, n;
35   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;
36 
37   PetscFunctionBegin;
38   PetscCall(MatGetSize(A, &m, &n));
39   PetscCall(PetscArrayzero(reductions, n));
40   if (type == NORM_2) {
41     for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] += PetscAbsScalar(aij->a[i] * aij->a[i]);
42   } else if (type == NORM_1) {
43     for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] += PetscAbsScalar(aij->a[i]);
44   } else if (type == NORM_INFINITY) {
45     for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] = PetscMax(PetscAbsScalar(aij->a[i]), reductions[aij->j[i]]);
46   } else if (type == REDUCTION_SUM_REALPART || type == REDUCTION_MEAN_REALPART) {
47     for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] += PetscRealPart(aij->a[i]);
48   } else if (type == REDUCTION_SUM_IMAGINARYPART || type == REDUCTION_MEAN_IMAGINARYPART) {
49     for (i = 0; i < aij->i[m]; i++) reductions[aij->j[i]] += PetscImaginaryPart(aij->a[i]);
50   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Unknown reduction type");
51 
52   if (type == NORM_2) {
53     for (i = 0; i < n; i++) reductions[i] = PetscSqrtReal(reductions[i]);
54   } else if (type == REDUCTION_MEAN_REALPART || type == REDUCTION_MEAN_IMAGINARYPART) {
55     for (i = 0; i < n; i++) reductions[i] /= m;
56   }
57   PetscFunctionReturn(PETSC_SUCCESS);
58 }
59 
60 PetscErrorCode MatFindOffBlockDiagonalEntries_SeqAIJ(Mat A, IS *is)
61 {
62   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
63   PetscInt        i, m = A->rmap->n, cnt = 0, bs = A->rmap->bs;
64   const PetscInt *jj = a->j, *ii = a->i;
65   PetscInt       *rows;
66 
67   PetscFunctionBegin;
68   for (i = 0; i < m; i++) {
69     if ((ii[i] != ii[i + 1]) && ((jj[ii[i]] < bs * (i / bs)) || (jj[ii[i + 1] - 1] > bs * ((i + bs) / bs) - 1))) cnt++;
70   }
71   PetscCall(PetscMalloc1(cnt, &rows));
72   cnt = 0;
73   for (i = 0; i < m; i++) {
74     if ((ii[i] != ii[i + 1]) && ((jj[ii[i]] < bs * (i / bs)) || (jj[ii[i + 1] - 1] > bs * ((i + bs) / bs) - 1))) {
75       rows[cnt] = i;
76       cnt++;
77     }
78   }
79   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cnt, rows, PETSC_OWN_POINTER, is));
80   PetscFunctionReturn(PETSC_SUCCESS);
81 }
82 
83 PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat A, PetscInt *nrows, PetscInt **zrows)
84 {
85   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
86   const MatScalar *aa;
87   PetscInt         i, m = A->rmap->n, cnt = 0;
88   const PetscInt  *ii = a->i, *jj = a->j, *diag;
89   PetscInt        *rows;
90 
91   PetscFunctionBegin;
92   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
93   PetscCall(MatMarkDiagonal_SeqAIJ(A));
94   diag = a->diag;
95   for (i = 0; i < m; i++) {
96     if ((diag[i] >= ii[i + 1]) || (jj[diag[i]] != i) || (aa[diag[i]] == 0.0)) cnt++;
97   }
98   PetscCall(PetscMalloc1(cnt, &rows));
99   cnt = 0;
100   for (i = 0; i < m; i++) {
101     if ((diag[i] >= ii[i + 1]) || (jj[diag[i]] != i) || (aa[diag[i]] == 0.0)) rows[cnt++] = i;
102   }
103   *nrows = cnt;
104   *zrows = rows;
105   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
106   PetscFunctionReturn(PETSC_SUCCESS);
107 }
108 
109 PetscErrorCode MatFindZeroDiagonals_SeqAIJ(Mat A, IS *zrows)
110 {
111   PetscInt nrows, *rows;
112 
113   PetscFunctionBegin;
114   *zrows = NULL;
115   PetscCall(MatFindZeroDiagonals_SeqAIJ_Private(A, &nrows, &rows));
116   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), nrows, rows, PETSC_OWN_POINTER, zrows));
117   PetscFunctionReturn(PETSC_SUCCESS);
118 }
119 
120 PetscErrorCode MatFindNonzeroRows_SeqAIJ(Mat A, IS *keptrows)
121 {
122   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
123   const MatScalar *aa;
124   PetscInt         m = A->rmap->n, cnt = 0;
125   const PetscInt  *ii;
126   PetscInt         n, i, j, *rows;
127 
128   PetscFunctionBegin;
129   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
130   *keptrows = NULL;
131   ii        = a->i;
132   for (i = 0; i < m; i++) {
133     n = ii[i + 1] - ii[i];
134     if (!n) {
135       cnt++;
136       goto ok1;
137     }
138     for (j = ii[i]; j < ii[i + 1]; j++) {
139       if (aa[j] != 0.0) goto ok1;
140     }
141     cnt++;
142   ok1:;
143   }
144   if (!cnt) {
145     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
146     PetscFunctionReturn(PETSC_SUCCESS);
147   }
148   PetscCall(PetscMalloc1(A->rmap->n - cnt, &rows));
149   cnt = 0;
150   for (i = 0; i < m; i++) {
151     n = ii[i + 1] - ii[i];
152     if (!n) continue;
153     for (j = ii[i]; j < ii[i + 1]; j++) {
154       if (aa[j] != 0.0) {
155         rows[cnt++] = i;
156         break;
157       }
158     }
159   }
160   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
161   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cnt, rows, PETSC_OWN_POINTER, keptrows));
162   PetscFunctionReturn(PETSC_SUCCESS);
163 }
164 
165 PetscErrorCode MatDiagonalSet_SeqAIJ(Mat Y, Vec D, InsertMode is)
166 {
167   Mat_SeqAIJ        *aij = (Mat_SeqAIJ *)Y->data;
168   PetscInt           i, m = Y->rmap->n;
169   const PetscInt    *diag;
170   MatScalar         *aa;
171   const PetscScalar *v;
172   PetscBool          missing;
173 
174   PetscFunctionBegin;
175   if (Y->assembled) {
176     PetscCall(MatMissingDiagonal_SeqAIJ(Y, &missing, NULL));
177     if (!missing) {
178       diag = aij->diag;
179       PetscCall(VecGetArrayRead(D, &v));
180       PetscCall(MatSeqAIJGetArray(Y, &aa));
181       if (is == INSERT_VALUES) {
182         for (i = 0; i < m; i++) aa[diag[i]] = v[i];
183       } else {
184         for (i = 0; i < m; i++) aa[diag[i]] += v[i];
185       }
186       PetscCall(MatSeqAIJRestoreArray(Y, &aa));
187       PetscCall(VecRestoreArrayRead(D, &v));
188       PetscFunctionReturn(PETSC_SUCCESS);
189     }
190     PetscCall(MatSeqAIJInvalidateDiagonal(Y));
191   }
192   PetscCall(MatDiagonalSet_Default(Y, D, is));
193   PetscFunctionReturn(PETSC_SUCCESS);
194 }
195 
196 PetscErrorCode MatGetRowIJ_SeqAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *m, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
197 {
198   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
199   PetscInt    i, ishift;
200 
201   PetscFunctionBegin;
202   if (m) *m = A->rmap->n;
203   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);
204   ishift = 0;
205   if (symmetric && A->structurally_symmetric != PETSC_BOOL3_TRUE) {
206     PetscCall(MatToSymmetricIJ_SeqAIJ(A->rmap->n, a->i, a->j, PETSC_TRUE, ishift, oshift, (PetscInt **)ia, (PetscInt **)ja));
207   } else if (oshift == 1) {
208     PetscInt *tia;
209     PetscInt  nz = a->i[A->rmap->n];
210     /* malloc space and  add 1 to i and j indices */
211     PetscCall(PetscMalloc1(A->rmap->n + 1, &tia));
212     for (i = 0; i < A->rmap->n + 1; i++) tia[i] = a->i[i] + 1;
213     *ia = tia;
214     if (ja) {
215       PetscInt *tja;
216       PetscCall(PetscMalloc1(nz + 1, &tja));
217       for (i = 0; i < nz; i++) tja[i] = a->j[i] + 1;
218       *ja = tja;
219     }
220   } else {
221     *ia = a->i;
222     if (ja) *ja = a->j;
223   }
224   PetscFunctionReturn(PETSC_SUCCESS);
225 }
226 
227 PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
228 {
229   PetscFunctionBegin;
230   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);
231   if ((symmetric && A->structurally_symmetric != PETSC_BOOL3_TRUE) || oshift == 1) {
232     PetscCall(PetscFree(*ia));
233     if (ja) PetscCall(PetscFree(*ja));
234   }
235   PetscFunctionReturn(PETSC_SUCCESS);
236 }
237 
238 PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *nn, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
239 {
240   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
241   PetscInt    i, *collengths, *cia, *cja, n = A->cmap->n, m = A->rmap->n;
242   PetscInt    nz = a->i[m], row, *jj, mr, col;
243 
244   PetscFunctionBegin;
245   *nn = n;
246   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);
247   if (symmetric) {
248     PetscCall(MatToSymmetricIJ_SeqAIJ(A->rmap->n, a->i, a->j, PETSC_TRUE, 0, oshift, (PetscInt **)ia, (PetscInt **)ja));
249   } else {
250     PetscCall(PetscCalloc1(n, &collengths));
251     PetscCall(PetscMalloc1(n + 1, &cia));
252     PetscCall(PetscMalloc1(nz, &cja));
253     jj = a->j;
254     for (i = 0; i < nz; i++) collengths[jj[i]]++;
255     cia[0] = oshift;
256     for (i = 0; i < n; i++) cia[i + 1] = cia[i] + collengths[i];
257     PetscCall(PetscArrayzero(collengths, n));
258     jj = a->j;
259     for (row = 0; row < m; row++) {
260       mr = a->i[row + 1] - a->i[row];
261       for (i = 0; i < mr; i++) {
262         col = *jj++;
263 
264         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
265       }
266     }
267     PetscCall(PetscFree(collengths));
268     *ia = cia;
269     *ja = cja;
270   }
271   PetscFunctionReturn(PETSC_SUCCESS);
272 }
273 
274 PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
275 {
276   PetscFunctionBegin;
277   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);
278 
279   PetscCall(PetscFree(*ia));
280   PetscCall(PetscFree(*ja));
281   PetscFunctionReturn(PETSC_SUCCESS);
282 }
283 
284 /*
285  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
286  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
287  spidx[], index of a->a, to be used in MatTransposeColoringCreate_SeqAIJ() and MatFDColoringCreate_SeqXAIJ()
288 */
289 PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *nn, const PetscInt *ia[], const PetscInt *ja[], PetscInt *spidx[], PetscBool *done)
290 {
291   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
292   PetscInt        i, *collengths, *cia, *cja, n = A->cmap->n, m = A->rmap->n;
293   PetscInt        nz = a->i[m], row, mr, col, tmp;
294   PetscInt       *cspidx;
295   const PetscInt *jj;
296 
297   PetscFunctionBegin;
298   *nn = n;
299   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);
300 
301   PetscCall(PetscCalloc1(n, &collengths));
302   PetscCall(PetscMalloc1(n + 1, &cia));
303   PetscCall(PetscMalloc1(nz, &cja));
304   PetscCall(PetscMalloc1(nz, &cspidx));
305   jj = a->j;
306   for (i = 0; i < nz; i++) collengths[jj[i]]++;
307   cia[0] = oshift;
308   for (i = 0; i < n; i++) cia[i + 1] = cia[i] + collengths[i];
309   PetscCall(PetscArrayzero(collengths, n));
310   jj = a->j;
311   for (row = 0; row < m; row++) {
312     mr = a->i[row + 1] - a->i[row];
313     for (i = 0; i < mr; i++) {
314       col         = *jj++;
315       tmp         = cia[col] + collengths[col]++ - oshift;
316       cspidx[tmp] = a->i[row] + i; /* index of a->j */
317       cja[tmp]    = row + oshift;
318     }
319   }
320   PetscCall(PetscFree(collengths));
321   *ia    = cia;
322   *ja    = cja;
323   *spidx = cspidx;
324   PetscFunctionReturn(PETSC_SUCCESS);
325 }
326 
327 PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscInt *spidx[], PetscBool *done)
328 {
329   PetscFunctionBegin;
330   PetscCall(MatRestoreColumnIJ_SeqAIJ(A, oshift, symmetric, inodecompressed, n, ia, ja, done));
331   PetscCall(PetscFree(*spidx));
332   PetscFunctionReturn(PETSC_SUCCESS);
333 }
334 
335 PetscErrorCode MatSetValuesRow_SeqAIJ(Mat A, PetscInt row, const PetscScalar v[])
336 {
337   Mat_SeqAIJ  *a  = (Mat_SeqAIJ *)A->data;
338   PetscInt    *ai = a->i;
339   PetscScalar *aa;
340 
341   PetscFunctionBegin;
342   PetscCall(MatSeqAIJGetArray(A, &aa));
343   PetscCall(PetscArraycpy(aa + ai[row], v, ai[row + 1] - ai[row]));
344   PetscCall(MatSeqAIJRestoreArray(A, &aa));
345   PetscFunctionReturn(PETSC_SUCCESS);
346 }
347 
348 /*
349     MatSeqAIJSetValuesLocalFast - An optimized version of MatSetValuesLocal() for SeqAIJ matrices with several assumptions
350 
351       -   a single row of values is set with each call
352       -   no row or column indices are negative or (in error) larger than the number of rows or columns
353       -   the values are always added to the matrix, not set
354       -   no new locations are introduced in the nonzero structure of the matrix
355 
356      This does NOT assume the global column indices are sorted
357 
358 */
359 
360 #include <petsc/private/isimpl.h>
361 PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
362 {
363   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
364   PetscInt        low, high, t, row, nrow, i, col, l;
365   const PetscInt *rp, *ai = a->i, *ailen = a->ilen, *aj = a->j;
366   PetscInt        lastcol = -1;
367   MatScalar      *ap, value, *aa;
368   const PetscInt *ridx = A->rmap->mapping->indices, *cidx = A->cmap->mapping->indices;
369 
370   PetscFunctionBegin;
371   PetscCall(MatSeqAIJGetArray(A, &aa));
372   row  = ridx[im[0]];
373   rp   = aj + ai[row];
374   ap   = aa + ai[row];
375   nrow = ailen[row];
376   low  = 0;
377   high = nrow;
378   for (l = 0; l < n; l++) { /* loop over added columns */
379     col   = cidx[in[l]];
380     value = v[l];
381 
382     if (col <= lastcol) low = 0;
383     else high = nrow;
384     lastcol = col;
385     while (high - low > 5) {
386       t = (low + high) / 2;
387       if (rp[t] > col) high = t;
388       else low = t;
389     }
390     for (i = low; i < high; i++) {
391       if (rp[i] == col) {
392         ap[i] += value;
393         low = i + 1;
394         break;
395       }
396     }
397   }
398   PetscCall(MatSeqAIJRestoreArray(A, &aa));
399   return PETSC_SUCCESS;
400 }
401 
402 PetscErrorCode MatSetValues_SeqAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
403 {
404   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
405   PetscInt   *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N;
406   PetscInt   *imax = a->imax, *ai = a->i, *ailen = a->ilen;
407   PetscInt   *aj = a->j, nonew = a->nonew, lastcol = -1;
408   MatScalar  *ap = NULL, value = 0.0, *aa;
409   PetscBool   ignorezeroentries = a->ignorezeroentries;
410   PetscBool   roworiented       = a->roworiented;
411 
412   PetscFunctionBegin;
413   PetscCall(MatSeqAIJGetArray(A, &aa));
414   for (k = 0; k < m; k++) { /* loop over added rows */
415     row = im[k];
416     if (row < 0) continue;
417     PetscCheck(row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, A->rmap->n - 1);
418     rp = aj + ai[row];
419     if (!A->structure_only) ap = aa + ai[row];
420     rmax = imax[row];
421     nrow = ailen[row];
422     low  = 0;
423     high = nrow;
424     for (l = 0; l < n; l++) { /* loop over added columns */
425       if (in[l] < 0) continue;
426       PetscCheck(in[l] < A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, in[l], A->cmap->n - 1);
427       col = in[l];
428       if (v && !A->structure_only) value = roworiented ? v[l + k * n] : v[k + l * m];
429       if (!A->structure_only && value == 0.0 && ignorezeroentries && is == ADD_VALUES && row != col) continue;
430 
431       if (col <= lastcol) low = 0;
432       else high = nrow;
433       lastcol = col;
434       while (high - low > 5) {
435         t = (low + high) / 2;
436         if (rp[t] > col) high = t;
437         else low = t;
438       }
439       for (i = low; i < high; i++) {
440         if (rp[i] > col) break;
441         if (rp[i] == col) {
442           if (!A->structure_only) {
443             if (is == ADD_VALUES) {
444               ap[i] += value;
445               (void)PetscLogFlops(1.0);
446             } else ap[i] = value;
447           }
448           low = i + 1;
449           goto noinsert;
450         }
451       }
452       if (value == 0.0 && ignorezeroentries && row != col) goto noinsert;
453       if (nonew == 1) goto noinsert;
454       PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") in the matrix", row, col);
455       if (A->structure_only) {
456         MatSeqXAIJReallocateAIJ_structure_only(A, A->rmap->n, 1, nrow, row, col, rmax, ai, aj, rp, imax, nonew, MatScalar);
457       } else {
458         MatSeqXAIJReallocateAIJ(A, A->rmap->n, 1, nrow, row, col, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar);
459       }
460       N = nrow++ - 1;
461       a->nz++;
462       high++;
463       /* shift up all the later entries in this row */
464       PetscCall(PetscArraymove(rp + i + 1, rp + i, N - i + 1));
465       rp[i] = col;
466       if (!A->structure_only) {
467         PetscCall(PetscArraymove(ap + i + 1, ap + i, N - i + 1));
468         ap[i] = value;
469       }
470       low = i + 1;
471       A->nonzerostate++;
472     noinsert:;
473     }
474     ailen[row] = nrow;
475   }
476   PetscCall(MatSeqAIJRestoreArray(A, &aa));
477   PetscFunctionReturn(PETSC_SUCCESS);
478 }
479 
480 PetscErrorCode MatSetValues_SeqAIJ_SortedFullNoPreallocation(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
481 {
482   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
483   PetscInt   *rp, k, row;
484   PetscInt   *ai = a->i;
485   PetscInt   *aj = a->j;
486   MatScalar  *aa, *ap;
487 
488   PetscFunctionBegin;
489   PetscCheck(!A->was_assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot call on assembled matrix.");
490   PetscCheck(m * n + a->nz <= a->maxnz, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of entries in matrix will be larger than maximum nonzeros allocated for %" PetscInt_FMT " in MatSeqAIJSetTotalPreallocation()", a->maxnz);
491 
492   PetscCall(MatSeqAIJGetArray(A, &aa));
493   for (k = 0; k < m; k++) { /* loop over added rows */
494     row = im[k];
495     rp  = aj + ai[row];
496     ap  = aa + ai[row];
497 
498     PetscCall(PetscMemcpy(rp, in, n * sizeof(PetscInt)));
499     if (!A->structure_only) {
500       if (v) {
501         PetscCall(PetscMemcpy(ap, v, n * sizeof(PetscScalar)));
502         v += n;
503       } else {
504         PetscCall(PetscMemzero(ap, n * sizeof(PetscScalar)));
505       }
506     }
507     a->ilen[row]  = n;
508     a->imax[row]  = n;
509     a->i[row + 1] = a->i[row] + n;
510     a->nz += n;
511   }
512   PetscCall(MatSeqAIJRestoreArray(A, &aa));
513   PetscFunctionReturn(PETSC_SUCCESS);
514 }
515 
516 /*@
517     MatSeqAIJSetTotalPreallocation - Sets an upper bound on the total number of expected nonzeros in the matrix.
518 
519   Input Parameters:
520 +  A - the `MATSEQAIJ` matrix
521 -  nztotal - bound on the number of nonzeros
522 
523   Level: advanced
524 
525   Notes:
526     This can be called if you will be provided the matrix row by row (from row zero) with sorted column indices for each row.
527     Simply call `MatSetValues()` after this call to provide the matrix entries in the usual manner. This matrix may be used
528     as always with multiple matrix assemblies.
529 
530 .seealso: [](ch_matrices), `Mat`, `MatSetOption()`, `MAT_SORTED_FULL`, `MatSetValues()`, `MatSeqAIJSetPreallocation()`
531 @*/
532 
533 PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat A, PetscInt nztotal)
534 {
535   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
536 
537   PetscFunctionBegin;
538   PetscCall(PetscLayoutSetUp(A->rmap));
539   PetscCall(PetscLayoutSetUp(A->cmap));
540   a->maxnz = nztotal;
541   if (!a->imax) { PetscCall(PetscMalloc1(A->rmap->n, &a->imax)); }
542   if (!a->ilen) {
543     PetscCall(PetscMalloc1(A->rmap->n, &a->ilen));
544   } else {
545     PetscCall(PetscMemzero(a->ilen, A->rmap->n * sizeof(PetscInt)));
546   }
547 
548   /* allocate the matrix space */
549   if (A->structure_only) {
550     PetscCall(PetscMalloc1(nztotal, &a->j));
551     PetscCall(PetscMalloc1(A->rmap->n + 1, &a->i));
552   } else {
553     PetscCall(PetscMalloc3(nztotal, &a->a, nztotal, &a->j, A->rmap->n + 1, &a->i));
554   }
555   a->i[0] = 0;
556   if (A->structure_only) {
557     a->singlemalloc = PETSC_FALSE;
558     a->free_a       = PETSC_FALSE;
559   } else {
560     a->singlemalloc = PETSC_TRUE;
561     a->free_a       = PETSC_TRUE;
562   }
563   a->free_ij        = PETSC_TRUE;
564   A->ops->setvalues = MatSetValues_SeqAIJ_SortedFullNoPreallocation;
565   A->preallocated   = PETSC_TRUE;
566   PetscFunctionReturn(PETSC_SUCCESS);
567 }
568 
569 PetscErrorCode MatSetValues_SeqAIJ_SortedFull(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
570 {
571   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
572   PetscInt   *rp, k, row;
573   PetscInt   *ai = a->i, *ailen = a->ilen;
574   PetscInt   *aj = a->j;
575   MatScalar  *aa, *ap;
576 
577   PetscFunctionBegin;
578   PetscCall(MatSeqAIJGetArray(A, &aa));
579   for (k = 0; k < m; k++) { /* loop over added rows */
580     row = im[k];
581     PetscCheck(n <= a->imax[row], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Preallocation for row %" PetscInt_FMT " does not match number of columns provided", n);
582     rp = aj + ai[row];
583     ap = aa + ai[row];
584     if (!A->was_assembled) PetscCall(PetscMemcpy(rp, in, n * sizeof(PetscInt)));
585     if (!A->structure_only) {
586       if (v) {
587         PetscCall(PetscMemcpy(ap, v, n * sizeof(PetscScalar)));
588         v += n;
589       } else {
590         PetscCall(PetscMemzero(ap, n * sizeof(PetscScalar)));
591       }
592     }
593     ailen[row] = n;
594     a->nz += n;
595   }
596   PetscCall(MatSeqAIJRestoreArray(A, &aa));
597   PetscFunctionReturn(PETSC_SUCCESS);
598 }
599 
600 PetscErrorCode MatGetValues_SeqAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], PetscScalar v[])
601 {
602   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
603   PetscInt        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
604   PetscInt        *ai = a->i, *ailen = a->ilen;
605   const MatScalar *ap, *aa;
606 
607   PetscFunctionBegin;
608   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
609   for (k = 0; k < m; k++) { /* loop over rows */
610     row = im[k];
611     if (row < 0) {
612       v += n;
613       continue;
614     } /* negative row */
615     PetscCheck(row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, A->rmap->n - 1);
616     rp   = aj + ai[row];
617     ap   = aa + ai[row];
618     nrow = ailen[row];
619     for (l = 0; l < n; l++) { /* loop over columns */
620       if (in[l] < 0) {
621         v++;
622         continue;
623       } /* negative column */
624       PetscCheck(in[l] < A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, in[l], A->cmap->n - 1);
625       col  = in[l];
626       high = nrow;
627       low  = 0; /* assume unsorted */
628       while (high - low > 5) {
629         t = (low + high) / 2;
630         if (rp[t] > col) high = t;
631         else low = t;
632       }
633       for (i = low; i < high; i++) {
634         if (rp[i] > col) break;
635         if (rp[i] == col) {
636           *v++ = ap[i];
637           goto finished;
638         }
639       }
640       *v++ = 0.0;
641     finished:;
642     }
643   }
644   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
645   PetscFunctionReturn(PETSC_SUCCESS);
646 }
647 
648 PetscErrorCode MatView_SeqAIJ_Binary(Mat mat, PetscViewer viewer)
649 {
650   Mat_SeqAIJ        *A = (Mat_SeqAIJ *)mat->data;
651   const PetscScalar *av;
652   PetscInt           header[4], M, N, m, nz, i;
653   PetscInt          *rowlens;
654 
655   PetscFunctionBegin;
656   PetscCall(PetscViewerSetUp(viewer));
657 
658   M  = mat->rmap->N;
659   N  = mat->cmap->N;
660   m  = mat->rmap->n;
661   nz = A->nz;
662 
663   /* write matrix header */
664   header[0] = MAT_FILE_CLASSID;
665   header[1] = M;
666   header[2] = N;
667   header[3] = nz;
668   PetscCall(PetscViewerBinaryWrite(viewer, header, 4, PETSC_INT));
669 
670   /* fill in and store row lengths */
671   PetscCall(PetscMalloc1(m, &rowlens));
672   for (i = 0; i < m; i++) rowlens[i] = A->i[i + 1] - A->i[i];
673   PetscCall(PetscViewerBinaryWrite(viewer, rowlens, m, PETSC_INT));
674   PetscCall(PetscFree(rowlens));
675   /* store column indices */
676   PetscCall(PetscViewerBinaryWrite(viewer, A->j, nz, PETSC_INT));
677   /* store nonzero values */
678   PetscCall(MatSeqAIJGetArrayRead(mat, &av));
679   PetscCall(PetscViewerBinaryWrite(viewer, av, nz, PETSC_SCALAR));
680   PetscCall(MatSeqAIJRestoreArrayRead(mat, &av));
681 
682   /* write block size option to the viewer's .info file */
683   PetscCall(MatView_Binary_BlockSizes(mat, viewer));
684   PetscFunctionReturn(PETSC_SUCCESS);
685 }
686 
687 static PetscErrorCode MatView_SeqAIJ_ASCII_structonly(Mat A, PetscViewer viewer)
688 {
689   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
690   PetscInt    i, k, m = A->rmap->N;
691 
692   PetscFunctionBegin;
693   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
694   for (i = 0; i < m; i++) {
695     PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
696     for (k = a->i[i]; k < a->i[i + 1]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ") ", a->j[k]));
697     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
698   }
699   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
700   PetscFunctionReturn(PETSC_SUCCESS);
701 }
702 
703 extern PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat, PetscViewer);
704 
705 PetscErrorCode MatView_SeqAIJ_ASCII(Mat A, PetscViewer viewer)
706 {
707   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
708   const PetscScalar *av;
709   PetscInt           i, j, m = A->rmap->n;
710   const char        *name;
711   PetscViewerFormat  format;
712 
713   PetscFunctionBegin;
714   if (A->structure_only) {
715     PetscCall(MatView_SeqAIJ_ASCII_structonly(A, viewer));
716     PetscFunctionReturn(PETSC_SUCCESS);
717   }
718 
719   PetscCall(PetscViewerGetFormat(viewer, &format));
720   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS);
721 
722   /* trigger copy to CPU if needed */
723   PetscCall(MatSeqAIJGetArrayRead(A, &av));
724   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
725   if (format == PETSC_VIEWER_ASCII_MATLAB) {
726     PetscInt nofinalvalue = 0;
727     if (m && ((a->i[m] == a->i[m - 1]) || (a->j[a->nz - 1] != A->cmap->n - 1))) {
728       /* Need a dummy value to ensure the dimension of the matrix. */
729       nofinalvalue = 1;
730     }
731     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
732     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Size = %" PetscInt_FMT " %" PetscInt_FMT " \n", m, A->cmap->n));
733     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Nonzeros = %" PetscInt_FMT " \n", a->nz));
734 #if defined(PETSC_USE_COMPLEX)
735     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",4);\n", a->nz + nofinalvalue));
736 #else
737     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",3);\n", a->nz + nofinalvalue));
738 #endif
739     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = [\n"));
740 
741     for (i = 0; i < m; i++) {
742       for (j = a->i[i]; j < a->i[i + 1]; j++) {
743 #if defined(PETSC_USE_COMPLEX)
744         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e %18.16e\n", i + 1, a->j[j] + 1, (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
745 #else
746         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n", i + 1, a->j[j] + 1, (double)a->a[j]));
747 #endif
748       }
749     }
750     if (nofinalvalue) {
751 #if defined(PETSC_USE_COMPLEX)
752       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e %18.16e\n", m, A->cmap->n, 0., 0.));
753 #else
754       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n", m, A->cmap->n, 0.0));
755 #endif
756     }
757     PetscCall(PetscObjectGetName((PetscObject)A, &name));
758     PetscCall(PetscViewerASCIIPrintf(viewer, "];\n %s = spconvert(zzz);\n", name));
759     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
760   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
761     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
762     for (i = 0; i < m; i++) {
763       PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
764       for (j = a->i[i]; j < a->i[i + 1]; j++) {
765 #if defined(PETSC_USE_COMPLEX)
766         if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
767           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
768         } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
769           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)-PetscImaginaryPart(a->a[j])));
770         } else if (PetscRealPart(a->a[j]) != 0.0) {
771           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j])));
772         }
773 #else
774         if (a->a[j] != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j]));
775 #endif
776       }
777       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
778     }
779     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
780   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
781     PetscInt nzd = 0, fshift = 1, *sptr;
782     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
783     PetscCall(PetscMalloc1(m + 1, &sptr));
784     for (i = 0; i < m; i++) {
785       sptr[i] = nzd + 1;
786       for (j = a->i[i]; j < a->i[i + 1]; j++) {
787         if (a->j[j] >= i) {
788 #if defined(PETSC_USE_COMPLEX)
789           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
790 #else
791           if (a->a[j] != 0.0) nzd++;
792 #endif
793         }
794       }
795     }
796     sptr[m] = nzd + 1;
797     PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT "\n\n", m, nzd));
798     for (i = 0; i < m + 1; i += 6) {
799       if (i + 4 < m) {
800         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2], sptr[i + 3], sptr[i + 4], sptr[i + 5]));
801       } else if (i + 3 < m) {
802         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2], sptr[i + 3], sptr[i + 4]));
803       } else if (i + 2 < m) {
804         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2], sptr[i + 3]));
805       } else if (i + 1 < m) {
806         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2]));
807       } else if (i < m) {
808         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1]));
809       } else {
810         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "\n", sptr[i]));
811       }
812     }
813     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
814     PetscCall(PetscFree(sptr));
815     for (i = 0; i < m; i++) {
816       for (j = a->i[i]; j < a->i[i + 1]; j++) {
817         if (a->j[j] >= i) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " ", a->j[j] + fshift));
818       }
819       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
820     }
821     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
822     for (i = 0; i < m; i++) {
823       for (j = a->i[i]; j < a->i[i + 1]; j++) {
824         if (a->j[j] >= i) {
825 #if defined(PETSC_USE_COMPLEX)
826           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " %18.16e %18.16e ", (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
827 #else
828           if (a->a[j] != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " %18.16e ", (double)a->a[j]));
829 #endif
830         }
831       }
832       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
833     }
834     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
835   } else if (format == PETSC_VIEWER_ASCII_DENSE) {
836     PetscInt    cnt = 0, jcnt;
837     PetscScalar value;
838 #if defined(PETSC_USE_COMPLEX)
839     PetscBool realonly = PETSC_TRUE;
840 
841     for (i = 0; i < a->i[m]; i++) {
842       if (PetscImaginaryPart(a->a[i]) != 0.0) {
843         realonly = PETSC_FALSE;
844         break;
845       }
846     }
847 #endif
848 
849     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
850     for (i = 0; i < m; i++) {
851       jcnt = 0;
852       for (j = 0; j < A->cmap->n; j++) {
853         if (jcnt < a->i[i + 1] - a->i[i] && j == a->j[cnt]) {
854           value = a->a[cnt++];
855           jcnt++;
856         } else {
857           value = 0.0;
858         }
859 #if defined(PETSC_USE_COMPLEX)
860         if (realonly) {
861           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)PetscRealPart(value)));
862         } else {
863           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e+%7.5e i ", (double)PetscRealPart(value), (double)PetscImaginaryPart(value)));
864         }
865 #else
866         PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)value));
867 #endif
868       }
869       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
870     }
871     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
872   } else if (format == PETSC_VIEWER_ASCII_MATRIXMARKET) {
873     PetscInt fshift = 1;
874     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
875 #if defined(PETSC_USE_COMPLEX)
876     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate complex general\n"));
877 #else
878     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate real general\n"));
879 #endif
880     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", m, A->cmap->n, a->nz));
881     for (i = 0; i < m; i++) {
882       for (j = a->i[i]; j < a->i[i + 1]; j++) {
883 #if defined(PETSC_USE_COMPLEX)
884         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %g %g\n", i + fshift, a->j[j] + fshift, (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
885 #else
886         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %g\n", i + fshift, a->j[j] + fshift, (double)a->a[j]));
887 #endif
888       }
889     }
890     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
891   } else {
892     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
893     if (A->factortype) {
894       for (i = 0; i < m; i++) {
895         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
896         /* L part */
897         for (j = a->i[i]; j < a->i[i + 1]; j++) {
898 #if defined(PETSC_USE_COMPLEX)
899           if (PetscImaginaryPart(a->a[j]) > 0.0) {
900             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
901           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
902             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)(-PetscImaginaryPart(a->a[j]))));
903           } else {
904             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j])));
905           }
906 #else
907           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j]));
908 #endif
909         }
910         /* diagonal */
911         j = a->diag[i];
912 #if defined(PETSC_USE_COMPLEX)
913         if (PetscImaginaryPart(a->a[j]) > 0.0) {
914           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(1.0 / a->a[j]), (double)PetscImaginaryPart(1.0 / a->a[j])));
915         } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
916           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(1.0 / a->a[j]), (double)(-PetscImaginaryPart(1.0 / a->a[j]))));
917         } else {
918           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(1.0 / a->a[j])));
919         }
920 #else
921         PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)(1.0 / a->a[j])));
922 #endif
923 
924         /* U part */
925         for (j = a->diag[i + 1] + 1; j < a->diag[i]; j++) {
926 #if defined(PETSC_USE_COMPLEX)
927           if (PetscImaginaryPart(a->a[j]) > 0.0) {
928             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
929           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
930             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)(-PetscImaginaryPart(a->a[j]))));
931           } else {
932             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j])));
933           }
934 #else
935           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j]));
936 #endif
937         }
938         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
939       }
940     } else {
941       for (i = 0; i < m; i++) {
942         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
943         for (j = a->i[i]; j < a->i[i + 1]; j++) {
944 #if defined(PETSC_USE_COMPLEX)
945           if (PetscImaginaryPart(a->a[j]) > 0.0) {
946             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
947           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
948             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)-PetscImaginaryPart(a->a[j])));
949           } else {
950             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j])));
951           }
952 #else
953           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j]));
954 #endif
955         }
956         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
957       }
958     }
959     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
960   }
961   PetscCall(PetscViewerFlush(viewer));
962   PetscFunctionReturn(PETSC_SUCCESS);
963 }
964 
965 #include <petscdraw.h>
966 PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw, void *Aa)
967 {
968   Mat                A = (Mat)Aa;
969   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
970   PetscInt           i, j, m = A->rmap->n;
971   int                color;
972   PetscReal          xl, yl, xr, yr, x_l, x_r, y_l, y_r;
973   PetscViewer        viewer;
974   PetscViewerFormat  format;
975   const PetscScalar *aa;
976 
977   PetscFunctionBegin;
978   PetscCall(PetscObjectQuery((PetscObject)A, "Zoomviewer", (PetscObject *)&viewer));
979   PetscCall(PetscViewerGetFormat(viewer, &format));
980   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
981 
982   /* loop over matrix elements drawing boxes */
983   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
984   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
985     PetscDrawCollectiveBegin(draw);
986     /* Blue for negative, Cyan for zero and  Red for positive */
987     color = PETSC_DRAW_BLUE;
988     for (i = 0; i < m; i++) {
989       y_l = m - i - 1.0;
990       y_r = y_l + 1.0;
991       for (j = a->i[i]; j < a->i[i + 1]; j++) {
992         x_l = a->j[j];
993         x_r = x_l + 1.0;
994         if (PetscRealPart(aa[j]) >= 0.) continue;
995         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
996       }
997     }
998     color = PETSC_DRAW_CYAN;
999     for (i = 0; i < m; i++) {
1000       y_l = m - i - 1.0;
1001       y_r = y_l + 1.0;
1002       for (j = a->i[i]; j < a->i[i + 1]; j++) {
1003         x_l = a->j[j];
1004         x_r = x_l + 1.0;
1005         if (aa[j] != 0.) continue;
1006         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1007       }
1008     }
1009     color = PETSC_DRAW_RED;
1010     for (i = 0; i < m; i++) {
1011       y_l = m - i - 1.0;
1012       y_r = y_l + 1.0;
1013       for (j = a->i[i]; j < a->i[i + 1]; j++) {
1014         x_l = a->j[j];
1015         x_r = x_l + 1.0;
1016         if (PetscRealPart(aa[j]) <= 0.) continue;
1017         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1018       }
1019     }
1020     PetscDrawCollectiveEnd(draw);
1021   } else {
1022     /* use contour shading to indicate magnitude of values */
1023     /* first determine max of all nonzero values */
1024     PetscReal minv = 0.0, maxv = 0.0;
1025     PetscInt  nz = a->nz, count = 0;
1026     PetscDraw popup;
1027 
1028     for (i = 0; i < nz; i++) {
1029       if (PetscAbsScalar(aa[i]) > maxv) maxv = PetscAbsScalar(aa[i]);
1030     }
1031     if (minv >= maxv) maxv = minv + PETSC_SMALL;
1032     PetscCall(PetscDrawGetPopup(draw, &popup));
1033     PetscCall(PetscDrawScalePopup(popup, minv, maxv));
1034 
1035     PetscDrawCollectiveBegin(draw);
1036     for (i = 0; i < m; i++) {
1037       y_l = m - i - 1.0;
1038       y_r = y_l + 1.0;
1039       for (j = a->i[i]; j < a->i[i + 1]; j++) {
1040         x_l   = a->j[j];
1041         x_r   = x_l + 1.0;
1042         color = PetscDrawRealToColor(PetscAbsScalar(aa[count]), minv, maxv);
1043         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1044         count++;
1045       }
1046     }
1047     PetscDrawCollectiveEnd(draw);
1048   }
1049   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1050   PetscFunctionReturn(PETSC_SUCCESS);
1051 }
1052 
1053 #include <petscdraw.h>
1054 PetscErrorCode MatView_SeqAIJ_Draw(Mat A, PetscViewer viewer)
1055 {
1056   PetscDraw draw;
1057   PetscReal xr, yr, xl, yl, h, w;
1058   PetscBool isnull;
1059 
1060   PetscFunctionBegin;
1061   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1062   PetscCall(PetscDrawIsNull(draw, &isnull));
1063   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
1064 
1065   xr = A->cmap->n;
1066   yr = A->rmap->n;
1067   h  = yr / 10.0;
1068   w  = xr / 10.0;
1069   xr += w;
1070   yr += h;
1071   xl = -w;
1072   yl = -h;
1073   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
1074   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", (PetscObject)viewer));
1075   PetscCall(PetscDrawZoom(draw, MatView_SeqAIJ_Draw_Zoom, A));
1076   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", NULL));
1077   PetscCall(PetscDrawSave(draw));
1078   PetscFunctionReturn(PETSC_SUCCESS);
1079 }
1080 
1081 PetscErrorCode MatView_SeqAIJ(Mat A, PetscViewer viewer)
1082 {
1083   PetscBool iascii, isbinary, isdraw;
1084 
1085   PetscFunctionBegin;
1086   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1087   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1088   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1089   if (iascii) PetscCall(MatView_SeqAIJ_ASCII(A, viewer));
1090   else if (isbinary) PetscCall(MatView_SeqAIJ_Binary(A, viewer));
1091   else if (isdraw) PetscCall(MatView_SeqAIJ_Draw(A, viewer));
1092   PetscCall(MatView_SeqAIJ_Inode(A, viewer));
1093   PetscFunctionReturn(PETSC_SUCCESS);
1094 }
1095 
1096 PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A, MatAssemblyType mode)
1097 {
1098   Mat_SeqAIJ *a      = (Mat_SeqAIJ *)A->data;
1099   PetscInt    fshift = 0, i, *ai = a->i, *aj = a->j, *imax = a->imax;
1100   PetscInt    m = A->rmap->n, *ip, N, *ailen = a->ilen, rmax = 0;
1101   MatScalar  *aa    = a->a, *ap;
1102   PetscReal   ratio = 0.6;
1103 
1104   PetscFunctionBegin;
1105   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
1106   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1107   if (A->was_assembled && A->ass_nonzerostate == A->nonzerostate) {
1108     /* we need to respect users asking to use or not the inodes routine in between matrix assemblies */
1109     PetscCall(MatAssemblyEnd_SeqAIJ_Inode(A, mode));
1110     PetscFunctionReturn(PETSC_SUCCESS);
1111   }
1112 
1113   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
1114   for (i = 1; i < m; i++) {
1115     /* move each row back by the amount of empty slots (fshift) before it*/
1116     fshift += imax[i - 1] - ailen[i - 1];
1117     rmax = PetscMax(rmax, ailen[i]);
1118     if (fshift) {
1119       ip = aj + ai[i];
1120       ap = aa + ai[i];
1121       N  = ailen[i];
1122       PetscCall(PetscArraymove(ip - fshift, ip, N));
1123       if (!A->structure_only) PetscCall(PetscArraymove(ap - fshift, ap, N));
1124     }
1125     ai[i] = ai[i - 1] + ailen[i - 1];
1126   }
1127   if (m) {
1128     fshift += imax[m - 1] - ailen[m - 1];
1129     ai[m] = ai[m - 1] + ailen[m - 1];
1130   }
1131   /* reset ilen and imax for each row */
1132   a->nonzerorowcnt = 0;
1133   if (A->structure_only) {
1134     PetscCall(PetscFree(a->imax));
1135     PetscCall(PetscFree(a->ilen));
1136   } else { /* !A->structure_only */
1137     for (i = 0; i < m; i++) {
1138       ailen[i] = imax[i] = ai[i + 1] - ai[i];
1139       a->nonzerorowcnt += ((ai[i + 1] - ai[i]) > 0);
1140     }
1141   }
1142   a->nz = ai[m];
1143   PetscCheck(!fshift || a->nounused != -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unused space detected in matrix: %" PetscInt_FMT " X %" PetscInt_FMT ", %" PetscInt_FMT " unneeded", m, A->cmap->n, fshift);
1144 
1145   PetscCall(MatMarkDiagonal_SeqAIJ(A));
1146   PetscCall(PetscInfo(A, "Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT "; storage space: %" PetscInt_FMT " unneeded,%" PetscInt_FMT " used\n", m, A->cmap->n, fshift, a->nz));
1147   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs));
1148   PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", rmax));
1149 
1150   A->info.mallocs += a->reallocs;
1151   a->reallocs         = 0;
1152   A->info.nz_unneeded = (PetscReal)fshift;
1153   a->rmax             = rmax;
1154 
1155   if (!A->structure_only) PetscCall(MatCheckCompressedRow(A, a->nonzerorowcnt, &a->compressedrow, a->i, m, ratio));
1156   PetscCall(MatAssemblyEnd_SeqAIJ_Inode(A, mode));
1157   PetscFunctionReturn(PETSC_SUCCESS);
1158 }
1159 
1160 PetscErrorCode MatRealPart_SeqAIJ(Mat A)
1161 {
1162   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1163   PetscInt    i, nz = a->nz;
1164   MatScalar  *aa;
1165 
1166   PetscFunctionBegin;
1167   PetscCall(MatSeqAIJGetArray(A, &aa));
1168   for (i = 0; i < nz; i++) aa[i] = PetscRealPart(aa[i]);
1169   PetscCall(MatSeqAIJRestoreArray(A, &aa));
1170   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1171   PetscFunctionReturn(PETSC_SUCCESS);
1172 }
1173 
1174 PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
1175 {
1176   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1177   PetscInt    i, nz = a->nz;
1178   MatScalar  *aa;
1179 
1180   PetscFunctionBegin;
1181   PetscCall(MatSeqAIJGetArray(A, &aa));
1182   for (i = 0; i < nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1183   PetscCall(MatSeqAIJRestoreArray(A, &aa));
1184   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1185   PetscFunctionReturn(PETSC_SUCCESS);
1186 }
1187 
1188 PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
1189 {
1190   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1191   MatScalar  *aa;
1192 
1193   PetscFunctionBegin;
1194   PetscCall(MatSeqAIJGetArrayWrite(A, &aa));
1195   PetscCall(PetscArrayzero(aa, a->i[A->rmap->n]));
1196   PetscCall(MatSeqAIJRestoreArrayWrite(A, &aa));
1197   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1198   PetscFunctionReturn(PETSC_SUCCESS);
1199 }
1200 
1201 PetscErrorCode MatDestroy_SeqAIJ(Mat A)
1202 {
1203   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1204 
1205   PetscFunctionBegin;
1206 #if defined(PETSC_USE_LOG)
1207   PetscCall(PetscLogObjectState((PetscObject)A, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, A->rmap->n, A->cmap->n, a->nz));
1208 #endif
1209   PetscCall(MatSeqXAIJFreeAIJ(A, &a->a, &a->j, &a->i));
1210   PetscCall(ISDestroy(&a->row));
1211   PetscCall(ISDestroy(&a->col));
1212   PetscCall(PetscFree(a->diag));
1213   PetscCall(PetscFree(a->ibdiag));
1214   PetscCall(PetscFree(a->imax));
1215   PetscCall(PetscFree(a->ilen));
1216   PetscCall(PetscFree(a->ipre));
1217   PetscCall(PetscFree3(a->idiag, a->mdiag, a->ssor_work));
1218   PetscCall(PetscFree(a->solve_work));
1219   PetscCall(ISDestroy(&a->icol));
1220   PetscCall(PetscFree(a->saved_values));
1221   PetscCall(PetscFree2(a->compressedrow.i, a->compressedrow.rindex));
1222   PetscCall(MatDestroy_SeqAIJ_Inode(A));
1223   PetscCall(PetscFree(A->data));
1224 
1225   /* MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted may allocate this.
1226      That function is so heavily used (sometimes in an hidden way through multnumeric function pointers)
1227      that is hard to properly add this data to the MatProduct data. We free it here to avoid
1228      users reusing the matrix object with different data to incur in obscure segmentation faults
1229      due to different matrix sizes */
1230   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc__ab_dense", NULL));
1231 
1232   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
1233   PetscCall(PetscObjectComposeFunction((PetscObject)A, "PetscMatlabEnginePut_C", NULL));
1234   PetscCall(PetscObjectComposeFunction((PetscObject)A, "PetscMatlabEngineGet_C", NULL));
1235   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJSetColumnIndices_C", NULL));
1236   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatStoreValues_C", NULL));
1237   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatRetrieveValues_C", NULL));
1238   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqsbaij_C", NULL));
1239   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqbaij_C", NULL));
1240   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijperm_C", NULL));
1241   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijsell_C", NULL));
1242 #if defined(PETSC_HAVE_MKL_SPARSE)
1243   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijmkl_C", NULL));
1244 #endif
1245 #if defined(PETSC_HAVE_CUDA)
1246   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijcusparse_C", NULL));
1247   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijcusparse_seqaij_C", NULL));
1248   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaij_seqaijcusparse_C", NULL));
1249 #endif
1250 #if defined(PETSC_HAVE_HIP)
1251   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijhipsparse_C", NULL));
1252   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijhipsparse_seqaij_C", NULL));
1253   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaij_seqaijhipsparse_C", NULL));
1254 #endif
1255 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1256   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijkokkos_C", NULL));
1257 #endif
1258   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijcrl_C", NULL));
1259 #if defined(PETSC_HAVE_ELEMENTAL)
1260   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_elemental_C", NULL));
1261 #endif
1262 #if defined(PETSC_HAVE_SCALAPACK)
1263   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_scalapack_C", NULL));
1264 #endif
1265 #if defined(PETSC_HAVE_HYPRE)
1266   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_hypre_C", NULL));
1267   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_transpose_seqaij_seqaij_C", NULL));
1268 #endif
1269   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqdense_C", NULL));
1270   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqsell_C", NULL));
1271   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_is_C", NULL));
1272   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatIsTranspose_C", NULL));
1273   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatIsHermitianTranspose_C", NULL));
1274   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", NULL));
1275   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatResetPreallocation_C", NULL));
1276   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJSetPreallocationCSR_C", NULL));
1277   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatReorderForNonzeroDiagonal_C", NULL));
1278   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_is_seqaij_C", NULL));
1279   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqdense_seqaij_C", NULL));
1280   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaij_seqaij_C", NULL));
1281   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJKron_C", NULL));
1282   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetPreallocationCOO_C", NULL));
1283   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetValuesCOO_C", NULL));
1284   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
1285   /* these calls do not belong here: the subclasses Duplicate/Destroy are wrong */
1286   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaijsell_seqaij_C", NULL));
1287   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaijperm_seqaij_C", NULL));
1288   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijviennacl_C", NULL));
1289   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijviennacl_seqdense_C", NULL));
1290   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijviennacl_seqaij_C", NULL));
1291   PetscFunctionReturn(PETSC_SUCCESS);
1292 }
1293 
1294 PetscErrorCode MatSetOption_SeqAIJ(Mat A, MatOption op, PetscBool flg)
1295 {
1296   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1297 
1298   PetscFunctionBegin;
1299   switch (op) {
1300   case MAT_ROW_ORIENTED:
1301     a->roworiented = flg;
1302     break;
1303   case MAT_KEEP_NONZERO_PATTERN:
1304     a->keepnonzeropattern = flg;
1305     break;
1306   case MAT_NEW_NONZERO_LOCATIONS:
1307     a->nonew = (flg ? 0 : 1);
1308     break;
1309   case MAT_NEW_NONZERO_LOCATION_ERR:
1310     a->nonew = (flg ? -1 : 0);
1311     break;
1312   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1313     a->nonew = (flg ? -2 : 0);
1314     break;
1315   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1316     a->nounused = (flg ? -1 : 0);
1317     break;
1318   case MAT_IGNORE_ZERO_ENTRIES:
1319     a->ignorezeroentries = flg;
1320     break;
1321   case MAT_SPD:
1322   case MAT_SYMMETRIC:
1323   case MAT_STRUCTURALLY_SYMMETRIC:
1324   case MAT_HERMITIAN:
1325   case MAT_SYMMETRY_ETERNAL:
1326   case MAT_STRUCTURE_ONLY:
1327   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
1328   case MAT_SPD_ETERNAL:
1329     /* if the diagonal matrix is square it inherits some of the properties above */
1330     break;
1331   case MAT_FORCE_DIAGONAL_ENTRIES:
1332   case MAT_IGNORE_OFF_PROC_ENTRIES:
1333   case MAT_USE_HASH_TABLE:
1334     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
1335     break;
1336   case MAT_USE_INODES:
1337     PetscCall(MatSetOption_SeqAIJ_Inode(A, MAT_USE_INODES, flg));
1338     break;
1339   case MAT_SUBMAT_SINGLEIS:
1340     A->submat_singleis = flg;
1341     break;
1342   case MAT_SORTED_FULL:
1343     if (flg) A->ops->setvalues = MatSetValues_SeqAIJ_SortedFull;
1344     else A->ops->setvalues = MatSetValues_SeqAIJ;
1345     break;
1346   case MAT_FORM_EXPLICIT_TRANSPOSE:
1347     A->form_explicit_transpose = flg;
1348     break;
1349   default:
1350     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
1351   }
1352   PetscFunctionReturn(PETSC_SUCCESS);
1353 }
1354 
1355 PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A, Vec v)
1356 {
1357   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1358   PetscInt           i, j, n, *ai = a->i, *aj = a->j;
1359   PetscScalar       *x;
1360   const PetscScalar *aa;
1361 
1362   PetscFunctionBegin;
1363   PetscCall(VecGetLocalSize(v, &n));
1364   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
1365   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
1366   if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) {
1367     PetscInt *diag = a->diag;
1368     PetscCall(VecGetArrayWrite(v, &x));
1369     for (i = 0; i < n; i++) x[i] = 1.0 / aa[diag[i]];
1370     PetscCall(VecRestoreArrayWrite(v, &x));
1371     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1372     PetscFunctionReturn(PETSC_SUCCESS);
1373   }
1374 
1375   PetscCall(VecGetArrayWrite(v, &x));
1376   for (i = 0; i < n; i++) {
1377     x[i] = 0.0;
1378     for (j = ai[i]; j < ai[i + 1]; j++) {
1379       if (aj[j] == i) {
1380         x[i] = aa[j];
1381         break;
1382       }
1383     }
1384   }
1385   PetscCall(VecRestoreArrayWrite(v, &x));
1386   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1387   PetscFunctionReturn(PETSC_SUCCESS);
1388 }
1389 
1390 #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>
1391 PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A, Vec xx, Vec zz, Vec yy)
1392 {
1393   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1394   const MatScalar   *aa;
1395   PetscScalar       *y;
1396   const PetscScalar *x;
1397   PetscInt           m = A->rmap->n;
1398 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1399   const MatScalar  *v;
1400   PetscScalar       alpha;
1401   PetscInt          n, i, j;
1402   const PetscInt   *idx, *ii, *ridx = NULL;
1403   Mat_CompressedRow cprow    = a->compressedrow;
1404   PetscBool         usecprow = cprow.use;
1405 #endif
1406 
1407   PetscFunctionBegin;
1408   if (zz != yy) PetscCall(VecCopy(zz, yy));
1409   PetscCall(VecGetArrayRead(xx, &x));
1410   PetscCall(VecGetArray(yy, &y));
1411   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
1412 
1413 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1414   fortranmulttransposeaddaij_(&m, x, a->i, a->j, aa, y);
1415 #else
1416   if (usecprow) {
1417     m    = cprow.nrows;
1418     ii   = cprow.i;
1419     ridx = cprow.rindex;
1420   } else {
1421     ii = a->i;
1422   }
1423   for (i = 0; i < m; i++) {
1424     idx = a->j + ii[i];
1425     v   = aa + ii[i];
1426     n   = ii[i + 1] - ii[i];
1427     if (usecprow) {
1428       alpha = x[ridx[i]];
1429     } else {
1430       alpha = x[i];
1431     }
1432     for (j = 0; j < n; j++) y[idx[j]] += alpha * v[j];
1433   }
1434 #endif
1435   PetscCall(PetscLogFlops(2.0 * a->nz));
1436   PetscCall(VecRestoreArrayRead(xx, &x));
1437   PetscCall(VecRestoreArray(yy, &y));
1438   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1439   PetscFunctionReturn(PETSC_SUCCESS);
1440 }
1441 
1442 PetscErrorCode MatMultTranspose_SeqAIJ(Mat A, Vec xx, Vec yy)
1443 {
1444   PetscFunctionBegin;
1445   PetscCall(VecSet(yy, 0.0));
1446   PetscCall(MatMultTransposeAdd_SeqAIJ(A, xx, yy, yy));
1447   PetscFunctionReturn(PETSC_SUCCESS);
1448 }
1449 
1450 #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>
1451 
1452 PetscErrorCode MatMult_SeqAIJ(Mat A, Vec xx, Vec yy)
1453 {
1454   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1455   PetscScalar       *y;
1456   const PetscScalar *x;
1457   const MatScalar   *aa, *a_a;
1458   PetscInt           m = A->rmap->n;
1459   const PetscInt    *aj, *ii, *ridx = NULL;
1460   PetscInt           n, i;
1461   PetscScalar        sum;
1462   PetscBool          usecprow = a->compressedrow.use;
1463 
1464 #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1465   #pragma disjoint(*x, *y, *aa)
1466 #endif
1467 
1468   PetscFunctionBegin;
1469   if (a->inode.use && a->inode.checked) {
1470     PetscCall(MatMult_SeqAIJ_Inode(A, xx, yy));
1471     PetscFunctionReturn(PETSC_SUCCESS);
1472   }
1473   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1474   PetscCall(VecGetArrayRead(xx, &x));
1475   PetscCall(VecGetArray(yy, &y));
1476   ii = a->i;
1477   if (usecprow) { /* use compressed row format */
1478     PetscCall(PetscArrayzero(y, m));
1479     m    = a->compressedrow.nrows;
1480     ii   = a->compressedrow.i;
1481     ridx = a->compressedrow.rindex;
1482     for (i = 0; i < m; i++) {
1483       n   = ii[i + 1] - ii[i];
1484       aj  = a->j + ii[i];
1485       aa  = a_a + ii[i];
1486       sum = 0.0;
1487       PetscSparseDensePlusDot(sum, x, aa, aj, n);
1488       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1489       y[*ridx++] = sum;
1490     }
1491   } else { /* do not use compressed row format */
1492 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1493     aj = a->j;
1494     aa = a_a;
1495     fortranmultaij_(&m, x, ii, aj, aa, y);
1496 #else
1497     for (i = 0; i < m; i++) {
1498       n   = ii[i + 1] - ii[i];
1499       aj  = a->j + ii[i];
1500       aa  = a_a + ii[i];
1501       sum = 0.0;
1502       PetscSparseDensePlusDot(sum, x, aa, aj, n);
1503       y[i] = sum;
1504     }
1505 #endif
1506   }
1507   PetscCall(PetscLogFlops(2.0 * a->nz - a->nonzerorowcnt));
1508   PetscCall(VecRestoreArrayRead(xx, &x));
1509   PetscCall(VecRestoreArray(yy, &y));
1510   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1511   PetscFunctionReturn(PETSC_SUCCESS);
1512 }
1513 
1514 PetscErrorCode MatMultMax_SeqAIJ(Mat A, Vec xx, Vec yy)
1515 {
1516   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1517   PetscScalar       *y;
1518   const PetscScalar *x;
1519   const MatScalar   *aa, *a_a;
1520   PetscInt           m = A->rmap->n;
1521   const PetscInt    *aj, *ii, *ridx   = NULL;
1522   PetscInt           n, i, nonzerorow = 0;
1523   PetscScalar        sum;
1524   PetscBool          usecprow = a->compressedrow.use;
1525 
1526 #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1527   #pragma disjoint(*x, *y, *aa)
1528 #endif
1529 
1530   PetscFunctionBegin;
1531   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1532   PetscCall(VecGetArrayRead(xx, &x));
1533   PetscCall(VecGetArray(yy, &y));
1534   if (usecprow) { /* use compressed row format */
1535     m    = a->compressedrow.nrows;
1536     ii   = a->compressedrow.i;
1537     ridx = a->compressedrow.rindex;
1538     for (i = 0; i < m; i++) {
1539       n   = ii[i + 1] - ii[i];
1540       aj  = a->j + ii[i];
1541       aa  = a_a + ii[i];
1542       sum = 0.0;
1543       nonzerorow += (n > 0);
1544       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1545       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1546       y[*ridx++] = sum;
1547     }
1548   } else { /* do not use compressed row format */
1549     ii = a->i;
1550     for (i = 0; i < m; i++) {
1551       n   = ii[i + 1] - ii[i];
1552       aj  = a->j + ii[i];
1553       aa  = a_a + ii[i];
1554       sum = 0.0;
1555       nonzerorow += (n > 0);
1556       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1557       y[i] = sum;
1558     }
1559   }
1560   PetscCall(PetscLogFlops(2.0 * a->nz - nonzerorow));
1561   PetscCall(VecRestoreArrayRead(xx, &x));
1562   PetscCall(VecRestoreArray(yy, &y));
1563   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1564   PetscFunctionReturn(PETSC_SUCCESS);
1565 }
1566 
1567 PetscErrorCode MatMultAddMax_SeqAIJ(Mat A, Vec xx, Vec yy, Vec zz)
1568 {
1569   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1570   PetscScalar       *y, *z;
1571   const PetscScalar *x;
1572   const MatScalar   *aa, *a_a;
1573   PetscInt           m = A->rmap->n, *aj, *ii;
1574   PetscInt           n, i, *ridx = NULL;
1575   PetscScalar        sum;
1576   PetscBool          usecprow = a->compressedrow.use;
1577 
1578   PetscFunctionBegin;
1579   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1580   PetscCall(VecGetArrayRead(xx, &x));
1581   PetscCall(VecGetArrayPair(yy, zz, &y, &z));
1582   if (usecprow) { /* use compressed row format */
1583     if (zz != yy) PetscCall(PetscArraycpy(z, y, m));
1584     m    = a->compressedrow.nrows;
1585     ii   = a->compressedrow.i;
1586     ridx = a->compressedrow.rindex;
1587     for (i = 0; i < m; i++) {
1588       n   = ii[i + 1] - ii[i];
1589       aj  = a->j + ii[i];
1590       aa  = a_a + ii[i];
1591       sum = y[*ridx];
1592       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1593       z[*ridx++] = sum;
1594     }
1595   } else { /* do not use compressed row format */
1596     ii = a->i;
1597     for (i = 0; i < m; i++) {
1598       n   = ii[i + 1] - ii[i];
1599       aj  = a->j + ii[i];
1600       aa  = a_a + ii[i];
1601       sum = y[i];
1602       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1603       z[i] = sum;
1604     }
1605   }
1606   PetscCall(PetscLogFlops(2.0 * a->nz));
1607   PetscCall(VecRestoreArrayRead(xx, &x));
1608   PetscCall(VecRestoreArrayPair(yy, zz, &y, &z));
1609   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1610   PetscFunctionReturn(PETSC_SUCCESS);
1611 }
1612 
1613 #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h>
1614 PetscErrorCode MatMultAdd_SeqAIJ(Mat A, Vec xx, Vec yy, Vec zz)
1615 {
1616   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1617   PetscScalar       *y, *z;
1618   const PetscScalar *x;
1619   const MatScalar   *aa, *a_a;
1620   const PetscInt    *aj, *ii, *ridx = NULL;
1621   PetscInt           m = A->rmap->n, n, i;
1622   PetscScalar        sum;
1623   PetscBool          usecprow = a->compressedrow.use;
1624 
1625   PetscFunctionBegin;
1626   if (a->inode.use && a->inode.checked) {
1627     PetscCall(MatMultAdd_SeqAIJ_Inode(A, xx, yy, zz));
1628     PetscFunctionReturn(PETSC_SUCCESS);
1629   }
1630   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1631   PetscCall(VecGetArrayRead(xx, &x));
1632   PetscCall(VecGetArrayPair(yy, zz, &y, &z));
1633   if (usecprow) { /* use compressed row format */
1634     if (zz != yy) PetscCall(PetscArraycpy(z, y, m));
1635     m    = a->compressedrow.nrows;
1636     ii   = a->compressedrow.i;
1637     ridx = a->compressedrow.rindex;
1638     for (i = 0; i < m; i++) {
1639       n   = ii[i + 1] - ii[i];
1640       aj  = a->j + ii[i];
1641       aa  = a_a + ii[i];
1642       sum = y[*ridx];
1643       PetscSparseDensePlusDot(sum, x, aa, aj, n);
1644       z[*ridx++] = sum;
1645     }
1646   } else { /* do not use compressed row format */
1647     ii = a->i;
1648 #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1649     aj = a->j;
1650     aa = a_a;
1651     fortranmultaddaij_(&m, x, ii, aj, aa, y, z);
1652 #else
1653     for (i = 0; i < m; i++) {
1654       n   = ii[i + 1] - ii[i];
1655       aj  = a->j + ii[i];
1656       aa  = a_a + ii[i];
1657       sum = y[i];
1658       PetscSparseDensePlusDot(sum, x, aa, aj, n);
1659       z[i] = sum;
1660     }
1661 #endif
1662   }
1663   PetscCall(PetscLogFlops(2.0 * a->nz));
1664   PetscCall(VecRestoreArrayRead(xx, &x));
1665   PetscCall(VecRestoreArrayPair(yy, zz, &y, &z));
1666   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1667   PetscFunctionReturn(PETSC_SUCCESS);
1668 }
1669 
1670 /*
1671      Adds diagonal pointers to sparse matrix structure.
1672 */
1673 PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1674 {
1675   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1676   PetscInt    i, j, m = A->rmap->n;
1677   PetscBool   alreadySet = PETSC_TRUE;
1678 
1679   PetscFunctionBegin;
1680   if (!a->diag) {
1681     PetscCall(PetscMalloc1(m, &a->diag));
1682     alreadySet = PETSC_FALSE;
1683   }
1684   for (i = 0; i < A->rmap->n; i++) {
1685     /* If A's diagonal is already correctly set, this fast track enables cheap and repeated MatMarkDiagonal_SeqAIJ() calls */
1686     if (alreadySet) {
1687       PetscInt pos = a->diag[i];
1688       if (pos >= a->i[i] && pos < a->i[i + 1] && a->j[pos] == i) continue;
1689     }
1690 
1691     a->diag[i] = a->i[i + 1];
1692     for (j = a->i[i]; j < a->i[i + 1]; j++) {
1693       if (a->j[j] == i) {
1694         a->diag[i] = j;
1695         break;
1696       }
1697     }
1698   }
1699   PetscFunctionReturn(PETSC_SUCCESS);
1700 }
1701 
1702 PetscErrorCode MatShift_SeqAIJ(Mat A, PetscScalar v)
1703 {
1704   Mat_SeqAIJ     *a    = (Mat_SeqAIJ *)A->data;
1705   const PetscInt *diag = (const PetscInt *)a->diag;
1706   const PetscInt *ii   = (const PetscInt *)a->i;
1707   PetscInt        i, *mdiag = NULL;
1708   PetscInt        cnt = 0; /* how many diagonals are missing */
1709 
1710   PetscFunctionBegin;
1711   if (!A->preallocated || !a->nz) {
1712     PetscCall(MatSeqAIJSetPreallocation(A, 1, NULL));
1713     PetscCall(MatShift_Basic(A, v));
1714     PetscFunctionReturn(PETSC_SUCCESS);
1715   }
1716 
1717   if (a->diagonaldense) {
1718     cnt = 0;
1719   } else {
1720     PetscCall(PetscCalloc1(A->rmap->n, &mdiag));
1721     for (i = 0; i < A->rmap->n; i++) {
1722       if (i < A->cmap->n && diag[i] >= ii[i + 1]) { /* 'out of range' rows never have diagonals */
1723         cnt++;
1724         mdiag[i] = 1;
1725       }
1726     }
1727   }
1728   if (!cnt) {
1729     PetscCall(MatShift_Basic(A, v));
1730   } else {
1731     PetscScalar *olda = a->a; /* preserve pointers to current matrix nonzeros structure and values */
1732     PetscInt    *oldj = a->j, *oldi = a->i;
1733     PetscBool    singlemalloc = a->singlemalloc, free_a = a->free_a, free_ij = a->free_ij;
1734 
1735     a->a = NULL;
1736     a->j = NULL;
1737     a->i = NULL;
1738     /* increase the values in imax for each row where a diagonal is being inserted then reallocate the matrix data structures */
1739     for (i = 0; i < PetscMin(A->rmap->n, A->cmap->n); i++) a->imax[i] += mdiag[i];
1740     PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(A, 0, a->imax));
1741 
1742     /* copy old values into new matrix data structure */
1743     for (i = 0; i < A->rmap->n; i++) {
1744       PetscCall(MatSetValues(A, 1, &i, a->imax[i] - mdiag[i], &oldj[oldi[i]], &olda[oldi[i]], ADD_VALUES));
1745       if (i < A->cmap->n) PetscCall(MatSetValue(A, i, i, v, ADD_VALUES));
1746     }
1747     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1748     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1749     if (singlemalloc) {
1750       PetscCall(PetscFree3(olda, oldj, oldi));
1751     } else {
1752       if (free_a) PetscCall(PetscFree(olda));
1753       if (free_ij) PetscCall(PetscFree(oldj));
1754       if (free_ij) PetscCall(PetscFree(oldi));
1755     }
1756   }
1757   PetscCall(PetscFree(mdiag));
1758   a->diagonaldense = PETSC_TRUE;
1759   PetscFunctionReturn(PETSC_SUCCESS);
1760 }
1761 
1762 /*
1763      Checks for missing diagonals
1764 */
1765 PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A, PetscBool *missing, PetscInt *d)
1766 {
1767   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1768   PetscInt   *diag, *ii = a->i, i;
1769 
1770   PetscFunctionBegin;
1771   *missing = PETSC_FALSE;
1772   if (A->rmap->n > 0 && !ii) {
1773     *missing = PETSC_TRUE;
1774     if (d) *d = 0;
1775     PetscCall(PetscInfo(A, "Matrix has no entries therefore is missing diagonal\n"));
1776   } else {
1777     PetscInt n;
1778     n    = PetscMin(A->rmap->n, A->cmap->n);
1779     diag = a->diag;
1780     for (i = 0; i < n; i++) {
1781       if (diag[i] >= ii[i + 1]) {
1782         *missing = PETSC_TRUE;
1783         if (d) *d = i;
1784         PetscCall(PetscInfo(A, "Matrix is missing diagonal number %" PetscInt_FMT "\n", i));
1785         break;
1786       }
1787     }
1788   }
1789   PetscFunctionReturn(PETSC_SUCCESS);
1790 }
1791 
1792 #include <petscblaslapack.h>
1793 #include <petsc/private/kernels/blockinvert.h>
1794 
1795 /*
1796     Note that values is allocated externally by the PC and then passed into this routine
1797 */
1798 PetscErrorCode MatInvertVariableBlockDiagonal_SeqAIJ(Mat A, PetscInt nblocks, const PetscInt *bsizes, PetscScalar *diag)
1799 {
1800   PetscInt        n = A->rmap->n, i, ncnt = 0, *indx, j, bsizemax = 0, *v_pivots;
1801   PetscBool       allowzeropivot, zeropivotdetected = PETSC_FALSE;
1802   const PetscReal shift = 0.0;
1803   PetscInt        ipvt[5];
1804   PetscCount      flops = 0;
1805   PetscScalar     work[25], *v_work;
1806 
1807   PetscFunctionBegin;
1808   allowzeropivot = PetscNot(A->erroriffailure);
1809   for (i = 0; i < nblocks; i++) ncnt += bsizes[i];
1810   PetscCheck(ncnt == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Total blocksizes %" PetscInt_FMT " doesn't match number matrix rows %" PetscInt_FMT, ncnt, n);
1811   for (i = 0; i < nblocks; i++) bsizemax = PetscMax(bsizemax, bsizes[i]);
1812   PetscCall(PetscMalloc1(bsizemax, &indx));
1813   if (bsizemax > 7) PetscCall(PetscMalloc2(bsizemax, &v_work, bsizemax, &v_pivots));
1814   ncnt = 0;
1815   for (i = 0; i < nblocks; i++) {
1816     for (j = 0; j < bsizes[i]; j++) indx[j] = ncnt + j;
1817     PetscCall(MatGetValues(A, bsizes[i], indx, bsizes[i], indx, diag));
1818     switch (bsizes[i]) {
1819     case 1:
1820       *diag = 1.0 / (*diag);
1821       break;
1822     case 2:
1823       PetscCall(PetscKernel_A_gets_inverse_A_2(diag, shift, allowzeropivot, &zeropivotdetected));
1824       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1825       PetscCall(PetscKernel_A_gets_transpose_A_2(diag));
1826       break;
1827     case 3:
1828       PetscCall(PetscKernel_A_gets_inverse_A_3(diag, shift, allowzeropivot, &zeropivotdetected));
1829       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1830       PetscCall(PetscKernel_A_gets_transpose_A_3(diag));
1831       break;
1832     case 4:
1833       PetscCall(PetscKernel_A_gets_inverse_A_4(diag, shift, allowzeropivot, &zeropivotdetected));
1834       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1835       PetscCall(PetscKernel_A_gets_transpose_A_4(diag));
1836       break;
1837     case 5:
1838       PetscCall(PetscKernel_A_gets_inverse_A_5(diag, ipvt, work, shift, allowzeropivot, &zeropivotdetected));
1839       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1840       PetscCall(PetscKernel_A_gets_transpose_A_5(diag));
1841       break;
1842     case 6:
1843       PetscCall(PetscKernel_A_gets_inverse_A_6(diag, shift, allowzeropivot, &zeropivotdetected));
1844       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1845       PetscCall(PetscKernel_A_gets_transpose_A_6(diag));
1846       break;
1847     case 7:
1848       PetscCall(PetscKernel_A_gets_inverse_A_7(diag, shift, allowzeropivot, &zeropivotdetected));
1849       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1850       PetscCall(PetscKernel_A_gets_transpose_A_7(diag));
1851       break;
1852     default:
1853       PetscCall(PetscKernel_A_gets_inverse_A(bsizes[i], diag, v_pivots, v_work, allowzeropivot, &zeropivotdetected));
1854       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1855       PetscCall(PetscKernel_A_gets_transpose_A_N(diag, bsizes[i]));
1856     }
1857     ncnt += bsizes[i];
1858     diag += bsizes[i] * bsizes[i];
1859     flops += 2 * PetscPowInt(bsizes[i], 3) / 3;
1860   }
1861   PetscCall(PetscLogFlops(flops));
1862   if (bsizemax > 7) PetscCall(PetscFree2(v_work, v_pivots));
1863   PetscCall(PetscFree(indx));
1864   PetscFunctionReturn(PETSC_SUCCESS);
1865 }
1866 
1867 /*
1868    Negative shift indicates do not generate an error if there is a zero diagonal, just invert it anyways
1869 */
1870 PetscErrorCode MatInvertDiagonal_SeqAIJ(Mat A, PetscScalar omega, PetscScalar fshift)
1871 {
1872   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
1873   PetscInt         i, *diag, m = A->rmap->n;
1874   const MatScalar *v;
1875   PetscScalar     *idiag, *mdiag;
1876 
1877   PetscFunctionBegin;
1878   if (a->idiagvalid) PetscFunctionReturn(PETSC_SUCCESS);
1879   PetscCall(MatMarkDiagonal_SeqAIJ(A));
1880   diag = a->diag;
1881   if (!a->idiag) { PetscCall(PetscMalloc3(m, &a->idiag, m, &a->mdiag, m, &a->ssor_work)); }
1882 
1883   mdiag = a->mdiag;
1884   idiag = a->idiag;
1885   PetscCall(MatSeqAIJGetArrayRead(A, &v));
1886   if (omega == 1.0 && PetscRealPart(fshift) <= 0.0) {
1887     for (i = 0; i < m; i++) {
1888       mdiag[i] = v[diag[i]];
1889       if (!PetscAbsScalar(mdiag[i])) { /* zero diagonal */
1890         if (PetscRealPart(fshift)) {
1891           PetscCall(PetscInfo(A, "Zero diagonal on row %" PetscInt_FMT "\n", i));
1892           A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1893           A->factorerror_zeropivot_value = 0.0;
1894           A->factorerror_zeropivot_row   = i;
1895         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Zero diagonal on row %" PetscInt_FMT, i);
1896       }
1897       idiag[i] = 1.0 / v[diag[i]];
1898     }
1899     PetscCall(PetscLogFlops(m));
1900   } else {
1901     for (i = 0; i < m; i++) {
1902       mdiag[i] = v[diag[i]];
1903       idiag[i] = omega / (fshift + v[diag[i]]);
1904     }
1905     PetscCall(PetscLogFlops(2.0 * m));
1906   }
1907   a->idiagvalid = PETSC_TRUE;
1908   PetscCall(MatSeqAIJRestoreArrayRead(A, &v));
1909   PetscFunctionReturn(PETSC_SUCCESS);
1910 }
1911 
1912 #include <../src/mat/impls/aij/seq/ftn-kernels/frelax.h>
1913 PetscErrorCode MatSOR_SeqAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
1914 {
1915   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1916   PetscScalar       *x, d, sum, *t, scale;
1917   const MatScalar   *v, *idiag = NULL, *mdiag, *aa;
1918   const PetscScalar *b, *bs, *xb, *ts;
1919   PetscInt           n, m = A->rmap->n, i;
1920   const PetscInt    *idx, *diag;
1921 
1922   PetscFunctionBegin;
1923   if (a->inode.use && a->inode.checked && omega == 1.0 && fshift == 0.0) {
1924     PetscCall(MatSOR_SeqAIJ_Inode(A, bb, omega, flag, fshift, its, lits, xx));
1925     PetscFunctionReturn(PETSC_SUCCESS);
1926   }
1927   its = its * lits;
1928 
1929   if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
1930   if (!a->idiagvalid) PetscCall(MatInvertDiagonal_SeqAIJ(A, omega, fshift));
1931   a->fshift = fshift;
1932   a->omega  = omega;
1933 
1934   diag  = a->diag;
1935   t     = a->ssor_work;
1936   idiag = a->idiag;
1937   mdiag = a->mdiag;
1938 
1939   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
1940   PetscCall(VecGetArray(xx, &x));
1941   PetscCall(VecGetArrayRead(bb, &b));
1942   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1943   if (flag == SOR_APPLY_UPPER) {
1944     /* apply (U + D/omega) to the vector */
1945     bs = b;
1946     for (i = 0; i < m; i++) {
1947       d   = fshift + mdiag[i];
1948       n   = a->i[i + 1] - diag[i] - 1;
1949       idx = a->j + diag[i] + 1;
1950       v   = aa + diag[i] + 1;
1951       sum = b[i] * d / omega;
1952       PetscSparseDensePlusDot(sum, bs, v, idx, n);
1953       x[i] = sum;
1954     }
1955     PetscCall(VecRestoreArray(xx, &x));
1956     PetscCall(VecRestoreArrayRead(bb, &b));
1957     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1958     PetscCall(PetscLogFlops(a->nz));
1959     PetscFunctionReturn(PETSC_SUCCESS);
1960   }
1961 
1962   PetscCheck(flag != SOR_APPLY_LOWER, PETSC_COMM_SELF, PETSC_ERR_SUP, "SOR_APPLY_LOWER is not implemented");
1963   if (flag & SOR_EISENSTAT) {
1964     /* Let  A = L + U + D; where L is lower triangular,
1965     U is upper triangular, E = D/omega; This routine applies
1966 
1967             (L + E)^{-1} A (U + E)^{-1}
1968 
1969     to a vector efficiently using Eisenstat's trick.
1970     */
1971     scale = (2.0 / omega) - 1.0;
1972 
1973     /*  x = (E + U)^{-1} b */
1974     for (i = m - 1; i >= 0; i--) {
1975       n   = a->i[i + 1] - diag[i] - 1;
1976       idx = a->j + diag[i] + 1;
1977       v   = aa + diag[i] + 1;
1978       sum = b[i];
1979       PetscSparseDenseMinusDot(sum, x, v, idx, n);
1980       x[i] = sum * idiag[i];
1981     }
1982 
1983     /*  t = b - (2*E - D)x */
1984     v = aa;
1985     for (i = 0; i < m; i++) t[i] = b[i] - scale * (v[*diag++]) * x[i];
1986 
1987     /*  t = (E + L)^{-1}t */
1988     ts   = t;
1989     diag = a->diag;
1990     for (i = 0; i < m; i++) {
1991       n   = diag[i] - a->i[i];
1992       idx = a->j + a->i[i];
1993       v   = aa + a->i[i];
1994       sum = t[i];
1995       PetscSparseDenseMinusDot(sum, ts, v, idx, n);
1996       t[i] = sum * idiag[i];
1997       /*  x = x + t */
1998       x[i] += t[i];
1999     }
2000 
2001     PetscCall(PetscLogFlops(6.0 * m - 1 + 2.0 * a->nz));
2002     PetscCall(VecRestoreArray(xx, &x));
2003     PetscCall(VecRestoreArrayRead(bb, &b));
2004     PetscFunctionReturn(PETSC_SUCCESS);
2005   }
2006   if (flag & SOR_ZERO_INITIAL_GUESS) {
2007     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
2008       for (i = 0; i < m; i++) {
2009         n   = diag[i] - a->i[i];
2010         idx = a->j + a->i[i];
2011         v   = aa + a->i[i];
2012         sum = b[i];
2013         PetscSparseDenseMinusDot(sum, x, v, idx, n);
2014         t[i] = sum;
2015         x[i] = sum * idiag[i];
2016       }
2017       xb = t;
2018       PetscCall(PetscLogFlops(a->nz));
2019     } else xb = b;
2020     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
2021       for (i = m - 1; i >= 0; i--) {
2022         n   = a->i[i + 1] - diag[i] - 1;
2023         idx = a->j + diag[i] + 1;
2024         v   = aa + diag[i] + 1;
2025         sum = xb[i];
2026         PetscSparseDenseMinusDot(sum, x, v, idx, n);
2027         if (xb == b) {
2028           x[i] = sum * idiag[i];
2029         } else {
2030           x[i] = (1 - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
2031         }
2032       }
2033       PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
2034     }
2035     its--;
2036   }
2037   while (its--) {
2038     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
2039       for (i = 0; i < m; i++) {
2040         /* lower */
2041         n   = diag[i] - a->i[i];
2042         idx = a->j + a->i[i];
2043         v   = aa + a->i[i];
2044         sum = b[i];
2045         PetscSparseDenseMinusDot(sum, x, v, idx, n);
2046         t[i] = sum; /* save application of the lower-triangular part */
2047         /* upper */
2048         n   = a->i[i + 1] - diag[i] - 1;
2049         idx = a->j + diag[i] + 1;
2050         v   = aa + diag[i] + 1;
2051         PetscSparseDenseMinusDot(sum, x, v, idx, n);
2052         x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
2053       }
2054       xb = t;
2055       PetscCall(PetscLogFlops(2.0 * a->nz));
2056     } else xb = b;
2057     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
2058       for (i = m - 1; i >= 0; i--) {
2059         sum = xb[i];
2060         if (xb == b) {
2061           /* whole matrix (no checkpointing available) */
2062           n   = a->i[i + 1] - a->i[i];
2063           idx = a->j + a->i[i];
2064           v   = aa + a->i[i];
2065           PetscSparseDenseMinusDot(sum, x, v, idx, n);
2066           x[i] = (1. - omega) * x[i] + (sum + mdiag[i] * x[i]) * idiag[i];
2067         } else { /* lower-triangular part has been saved, so only apply upper-triangular */
2068           n   = a->i[i + 1] - diag[i] - 1;
2069           idx = a->j + diag[i] + 1;
2070           v   = aa + diag[i] + 1;
2071           PetscSparseDenseMinusDot(sum, x, v, idx, n);
2072           x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
2073         }
2074       }
2075       if (xb == b) {
2076         PetscCall(PetscLogFlops(2.0 * a->nz));
2077       } else {
2078         PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
2079       }
2080     }
2081   }
2082   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2083   PetscCall(VecRestoreArray(xx, &x));
2084   PetscCall(VecRestoreArrayRead(bb, &b));
2085   PetscFunctionReturn(PETSC_SUCCESS);
2086 }
2087 
2088 PetscErrorCode MatGetInfo_SeqAIJ(Mat A, MatInfoType flag, MatInfo *info)
2089 {
2090   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
2091 
2092   PetscFunctionBegin;
2093   info->block_size   = 1.0;
2094   info->nz_allocated = a->maxnz;
2095   info->nz_used      = a->nz;
2096   info->nz_unneeded  = (a->maxnz - a->nz);
2097   info->assemblies   = A->num_ass;
2098   info->mallocs      = A->info.mallocs;
2099   info->memory       = 0; /* REVIEW ME */
2100   if (A->factortype) {
2101     info->fill_ratio_given  = A->info.fill_ratio_given;
2102     info->fill_ratio_needed = A->info.fill_ratio_needed;
2103     info->factor_mallocs    = A->info.factor_mallocs;
2104   } else {
2105     info->fill_ratio_given  = 0;
2106     info->fill_ratio_needed = 0;
2107     info->factor_mallocs    = 0;
2108   }
2109   PetscFunctionReturn(PETSC_SUCCESS);
2110 }
2111 
2112 PetscErrorCode MatZeroRows_SeqAIJ(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
2113 {
2114   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2115   PetscInt           i, m = A->rmap->n - 1;
2116   const PetscScalar *xx;
2117   PetscScalar       *bb, *aa;
2118   PetscInt           d = 0;
2119 
2120   PetscFunctionBegin;
2121   if (x && b) {
2122     PetscCall(VecGetArrayRead(x, &xx));
2123     PetscCall(VecGetArray(b, &bb));
2124     for (i = 0; i < N; i++) {
2125       PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2126       if (rows[i] >= A->cmap->n) continue;
2127       bb[rows[i]] = diag * xx[rows[i]];
2128     }
2129     PetscCall(VecRestoreArrayRead(x, &xx));
2130     PetscCall(VecRestoreArray(b, &bb));
2131   }
2132 
2133   PetscCall(MatSeqAIJGetArray(A, &aa));
2134   if (a->keepnonzeropattern) {
2135     for (i = 0; i < N; i++) {
2136       PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2137       PetscCall(PetscArrayzero(&aa[a->i[rows[i]]], a->ilen[rows[i]]));
2138     }
2139     if (diag != 0.0) {
2140       for (i = 0; i < N; i++) {
2141         d = rows[i];
2142         if (rows[i] >= A->cmap->n) continue;
2143         PetscCheck(a->diag[d] < a->i[d + 1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry in the zeroed row %" PetscInt_FMT, d);
2144       }
2145       for (i = 0; i < N; i++) {
2146         if (rows[i] >= A->cmap->n) continue;
2147         aa[a->diag[rows[i]]] = diag;
2148       }
2149     }
2150   } else {
2151     if (diag != 0.0) {
2152       for (i = 0; i < N; i++) {
2153         PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2154         if (a->ilen[rows[i]] > 0) {
2155           if (rows[i] >= A->cmap->n) {
2156             a->ilen[rows[i]] = 0;
2157           } else {
2158             a->ilen[rows[i]]    = 1;
2159             aa[a->i[rows[i]]]   = diag;
2160             a->j[a->i[rows[i]]] = rows[i];
2161           }
2162         } else if (rows[i] < A->cmap->n) { /* in case row was completely empty */
2163           PetscCall(MatSetValues_SeqAIJ(A, 1, &rows[i], 1, &rows[i], &diag, INSERT_VALUES));
2164         }
2165       }
2166     } else {
2167       for (i = 0; i < N; i++) {
2168         PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2169         a->ilen[rows[i]] = 0;
2170       }
2171     }
2172     A->nonzerostate++;
2173   }
2174   PetscCall(MatSeqAIJRestoreArray(A, &aa));
2175   PetscUseTypeMethod(A, assemblyend, MAT_FINAL_ASSEMBLY);
2176   PetscFunctionReturn(PETSC_SUCCESS);
2177 }
2178 
2179 PetscErrorCode MatZeroRowsColumns_SeqAIJ(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
2180 {
2181   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2182   PetscInt           i, j, m = A->rmap->n - 1, d = 0;
2183   PetscBool          missing, *zeroed, vecs = PETSC_FALSE;
2184   const PetscScalar *xx;
2185   PetscScalar       *bb, *aa;
2186 
2187   PetscFunctionBegin;
2188   if (!N) PetscFunctionReturn(PETSC_SUCCESS);
2189   PetscCall(MatSeqAIJGetArray(A, &aa));
2190   if (x && b) {
2191     PetscCall(VecGetArrayRead(x, &xx));
2192     PetscCall(VecGetArray(b, &bb));
2193     vecs = PETSC_TRUE;
2194   }
2195   PetscCall(PetscCalloc1(A->rmap->n, &zeroed));
2196   for (i = 0; i < N; i++) {
2197     PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2198     PetscCall(PetscArrayzero(&aa[a->i[rows[i]]], a->ilen[rows[i]]));
2199 
2200     zeroed[rows[i]] = PETSC_TRUE;
2201   }
2202   for (i = 0; i < A->rmap->n; i++) {
2203     if (!zeroed[i]) {
2204       for (j = a->i[i]; j < a->i[i + 1]; j++) {
2205         if (a->j[j] < A->rmap->n && zeroed[a->j[j]]) {
2206           if (vecs) bb[i] -= aa[j] * xx[a->j[j]];
2207           aa[j] = 0.0;
2208         }
2209       }
2210     } else if (vecs && i < A->cmap->N) bb[i] = diag * xx[i];
2211   }
2212   if (x && b) {
2213     PetscCall(VecRestoreArrayRead(x, &xx));
2214     PetscCall(VecRestoreArray(b, &bb));
2215   }
2216   PetscCall(PetscFree(zeroed));
2217   if (diag != 0.0) {
2218     PetscCall(MatMissingDiagonal_SeqAIJ(A, &missing, &d));
2219     if (missing) {
2220       for (i = 0; i < N; i++) {
2221         if (rows[i] >= A->cmap->N) continue;
2222         PetscCheck(!a->nonew || rows[i] < d, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry in row %" PetscInt_FMT " (%" PetscInt_FMT ")", d, rows[i]);
2223         PetscCall(MatSetValues_SeqAIJ(A, 1, &rows[i], 1, &rows[i], &diag, INSERT_VALUES));
2224       }
2225     } else {
2226       for (i = 0; i < N; i++) aa[a->diag[rows[i]]] = diag;
2227     }
2228   }
2229   PetscCall(MatSeqAIJRestoreArray(A, &aa));
2230   PetscUseTypeMethod(A, assemblyend, MAT_FINAL_ASSEMBLY);
2231   PetscFunctionReturn(PETSC_SUCCESS);
2232 }
2233 
2234 PetscErrorCode MatGetRow_SeqAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
2235 {
2236   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2237   const PetscScalar *aa;
2238   PetscInt          *itmp;
2239 
2240   PetscFunctionBegin;
2241   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2242   *nz = a->i[row + 1] - a->i[row];
2243   if (v) *v = (PetscScalar *)(aa + a->i[row]);
2244   if (idx) {
2245     itmp = a->j + a->i[row];
2246     if (*nz) *idx = itmp;
2247     else *idx = NULL;
2248   }
2249   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2250   PetscFunctionReturn(PETSC_SUCCESS);
2251 }
2252 
2253 PetscErrorCode MatRestoreRow_SeqAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
2254 {
2255   PetscFunctionBegin;
2256   if (nz) *nz = 0;
2257   if (idx) *idx = NULL;
2258   if (v) *v = NULL;
2259   PetscFunctionReturn(PETSC_SUCCESS);
2260 }
2261 
2262 PetscErrorCode MatNorm_SeqAIJ(Mat A, NormType type, PetscReal *nrm)
2263 {
2264   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
2265   const MatScalar *v;
2266   PetscReal        sum = 0.0;
2267   PetscInt         i, j;
2268 
2269   PetscFunctionBegin;
2270   PetscCall(MatSeqAIJGetArrayRead(A, &v));
2271   if (type == NORM_FROBENIUS) {
2272 #if defined(PETSC_USE_REAL___FP16)
2273     PetscBLASInt one = 1, nz = a->nz;
2274     PetscCallBLAS("BLASnrm2", *nrm = BLASnrm2_(&nz, v, &one));
2275 #else
2276     for (i = 0; i < a->nz; i++) {
2277       sum += PetscRealPart(PetscConj(*v) * (*v));
2278       v++;
2279     }
2280     *nrm = PetscSqrtReal(sum);
2281 #endif
2282     PetscCall(PetscLogFlops(2.0 * a->nz));
2283   } else if (type == NORM_1) {
2284     PetscReal *tmp;
2285     PetscInt  *jj = a->j;
2286     PetscCall(PetscCalloc1(A->cmap->n + 1, &tmp));
2287     *nrm = 0.0;
2288     for (j = 0; j < a->nz; j++) {
2289       tmp[*jj++] += PetscAbsScalar(*v);
2290       v++;
2291     }
2292     for (j = 0; j < A->cmap->n; j++) {
2293       if (tmp[j] > *nrm) *nrm = tmp[j];
2294     }
2295     PetscCall(PetscFree(tmp));
2296     PetscCall(PetscLogFlops(PetscMax(a->nz - 1, 0)));
2297   } else if (type == NORM_INFINITY) {
2298     *nrm = 0.0;
2299     for (j = 0; j < A->rmap->n; j++) {
2300       const PetscScalar *v2 = v + a->i[j];
2301       sum                   = 0.0;
2302       for (i = 0; i < a->i[j + 1] - a->i[j]; i++) {
2303         sum += PetscAbsScalar(*v2);
2304         v2++;
2305       }
2306       if (sum > *nrm) *nrm = sum;
2307     }
2308     PetscCall(PetscLogFlops(PetscMax(a->nz - 1, 0)));
2309   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for two norm");
2310   PetscCall(MatSeqAIJRestoreArrayRead(A, &v));
2311   PetscFunctionReturn(PETSC_SUCCESS);
2312 }
2313 
2314 PetscErrorCode MatIsTranspose_SeqAIJ(Mat A, Mat B, PetscReal tol, PetscBool *f)
2315 {
2316   Mat_SeqAIJ      *aij = (Mat_SeqAIJ *)A->data, *bij = (Mat_SeqAIJ *)B->data;
2317   PetscInt        *adx, *bdx, *aii, *bii, *aptr, *bptr;
2318   const MatScalar *va, *vb;
2319   PetscInt         ma, na, mb, nb, i;
2320 
2321   PetscFunctionBegin;
2322   PetscCall(MatGetSize(A, &ma, &na));
2323   PetscCall(MatGetSize(B, &mb, &nb));
2324   if (ma != nb || na != mb) {
2325     *f = PETSC_FALSE;
2326     PetscFunctionReturn(PETSC_SUCCESS);
2327   }
2328   PetscCall(MatSeqAIJGetArrayRead(A, &va));
2329   PetscCall(MatSeqAIJGetArrayRead(B, &vb));
2330   aii = aij->i;
2331   bii = bij->i;
2332   adx = aij->j;
2333   bdx = bij->j;
2334   PetscCall(PetscMalloc1(ma, &aptr));
2335   PetscCall(PetscMalloc1(mb, &bptr));
2336   for (i = 0; i < ma; i++) aptr[i] = aii[i];
2337   for (i = 0; i < mb; i++) bptr[i] = bii[i];
2338 
2339   *f = PETSC_TRUE;
2340   for (i = 0; i < ma; i++) {
2341     while (aptr[i] < aii[i + 1]) {
2342       PetscInt    idc, idr;
2343       PetscScalar vc, vr;
2344       /* column/row index/value */
2345       idc = adx[aptr[i]];
2346       idr = bdx[bptr[idc]];
2347       vc  = va[aptr[i]];
2348       vr  = vb[bptr[idc]];
2349       if (i != idr || PetscAbsScalar(vc - vr) > tol) {
2350         *f = PETSC_FALSE;
2351         goto done;
2352       } else {
2353         aptr[i]++;
2354         if (B || i != idc) bptr[idc]++;
2355       }
2356     }
2357   }
2358 done:
2359   PetscCall(PetscFree(aptr));
2360   PetscCall(PetscFree(bptr));
2361   PetscCall(MatSeqAIJRestoreArrayRead(A, &va));
2362   PetscCall(MatSeqAIJRestoreArrayRead(B, &vb));
2363   PetscFunctionReturn(PETSC_SUCCESS);
2364 }
2365 
2366 PetscErrorCode MatIsHermitianTranspose_SeqAIJ(Mat A, Mat B, PetscReal tol, PetscBool *f)
2367 {
2368   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data, *bij = (Mat_SeqAIJ *)B->data;
2369   PetscInt   *adx, *bdx, *aii, *bii, *aptr, *bptr;
2370   MatScalar  *va, *vb;
2371   PetscInt    ma, na, mb, nb, i;
2372 
2373   PetscFunctionBegin;
2374   PetscCall(MatGetSize(A, &ma, &na));
2375   PetscCall(MatGetSize(B, &mb, &nb));
2376   if (ma != nb || na != mb) {
2377     *f = PETSC_FALSE;
2378     PetscFunctionReturn(PETSC_SUCCESS);
2379   }
2380   aii = aij->i;
2381   bii = bij->i;
2382   adx = aij->j;
2383   bdx = bij->j;
2384   va  = aij->a;
2385   vb  = bij->a;
2386   PetscCall(PetscMalloc1(ma, &aptr));
2387   PetscCall(PetscMalloc1(mb, &bptr));
2388   for (i = 0; i < ma; i++) aptr[i] = aii[i];
2389   for (i = 0; i < mb; i++) bptr[i] = bii[i];
2390 
2391   *f = PETSC_TRUE;
2392   for (i = 0; i < ma; i++) {
2393     while (aptr[i] < aii[i + 1]) {
2394       PetscInt    idc, idr;
2395       PetscScalar vc, vr;
2396       /* column/row index/value */
2397       idc = adx[aptr[i]];
2398       idr = bdx[bptr[idc]];
2399       vc  = va[aptr[i]];
2400       vr  = vb[bptr[idc]];
2401       if (i != idr || PetscAbsScalar(vc - PetscConj(vr)) > tol) {
2402         *f = PETSC_FALSE;
2403         goto done;
2404       } else {
2405         aptr[i]++;
2406         if (B || i != idc) bptr[idc]++;
2407       }
2408     }
2409   }
2410 done:
2411   PetscCall(PetscFree(aptr));
2412   PetscCall(PetscFree(bptr));
2413   PetscFunctionReturn(PETSC_SUCCESS);
2414 }
2415 
2416 PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A, PetscReal tol, PetscBool *f)
2417 {
2418   PetscFunctionBegin;
2419   PetscCall(MatIsTranspose_SeqAIJ(A, A, tol, f));
2420   PetscFunctionReturn(PETSC_SUCCESS);
2421 }
2422 
2423 PetscErrorCode MatIsHermitian_SeqAIJ(Mat A, PetscReal tol, PetscBool *f)
2424 {
2425   PetscFunctionBegin;
2426   PetscCall(MatIsHermitianTranspose_SeqAIJ(A, A, tol, f));
2427   PetscFunctionReturn(PETSC_SUCCESS);
2428 }
2429 
2430 PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A, Vec ll, Vec rr)
2431 {
2432   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2433   const PetscScalar *l, *r;
2434   PetscScalar        x;
2435   MatScalar         *v;
2436   PetscInt           i, j, m = A->rmap->n, n = A->cmap->n, M, nz = a->nz;
2437   const PetscInt    *jj;
2438 
2439   PetscFunctionBegin;
2440   if (ll) {
2441     /* The local size is used so that VecMPI can be passed to this routine
2442        by MatDiagonalScale_MPIAIJ */
2443     PetscCall(VecGetLocalSize(ll, &m));
2444     PetscCheck(m == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Left scaling vector wrong length");
2445     PetscCall(VecGetArrayRead(ll, &l));
2446     PetscCall(MatSeqAIJGetArray(A, &v));
2447     for (i = 0; i < m; i++) {
2448       x = l[i];
2449       M = a->i[i + 1] - a->i[i];
2450       for (j = 0; j < M; j++) (*v++) *= x;
2451     }
2452     PetscCall(VecRestoreArrayRead(ll, &l));
2453     PetscCall(PetscLogFlops(nz));
2454     PetscCall(MatSeqAIJRestoreArray(A, &v));
2455   }
2456   if (rr) {
2457     PetscCall(VecGetLocalSize(rr, &n));
2458     PetscCheck(n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Right scaling vector wrong length");
2459     PetscCall(VecGetArrayRead(rr, &r));
2460     PetscCall(MatSeqAIJGetArray(A, &v));
2461     jj = a->j;
2462     for (i = 0; i < nz; i++) (*v++) *= r[*jj++];
2463     PetscCall(MatSeqAIJRestoreArray(A, &v));
2464     PetscCall(VecRestoreArrayRead(rr, &r));
2465     PetscCall(PetscLogFlops(nz));
2466   }
2467   PetscCall(MatSeqAIJInvalidateDiagonal(A));
2468   PetscFunctionReturn(PETSC_SUCCESS);
2469 }
2470 
2471 PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat A, IS isrow, IS iscol, PetscInt csize, MatReuse scall, Mat *B)
2472 {
2473   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data, *c;
2474   PetscInt          *smap, i, k, kstart, kend, oldcols = A->cmap->n, *lens;
2475   PetscInt           row, mat_i, *mat_j, tcol, first, step, *mat_ilen, sum, lensi;
2476   const PetscInt    *irow, *icol;
2477   const PetscScalar *aa;
2478   PetscInt           nrows, ncols;
2479   PetscInt          *starts, *j_new, *i_new, *aj = a->j, *ai = a->i, ii, *ailen = a->ilen;
2480   MatScalar         *a_new, *mat_a, *c_a;
2481   Mat                C;
2482   PetscBool          stride;
2483 
2484   PetscFunctionBegin;
2485   PetscCall(ISGetIndices(isrow, &irow));
2486   PetscCall(ISGetLocalSize(isrow, &nrows));
2487   PetscCall(ISGetLocalSize(iscol, &ncols));
2488 
2489   PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISSTRIDE, &stride));
2490   if (stride) {
2491     PetscCall(ISStrideGetInfo(iscol, &first, &step));
2492   } else {
2493     first = 0;
2494     step  = 0;
2495   }
2496   if (stride && step == 1) {
2497     /* special case of contiguous rows */
2498     PetscCall(PetscMalloc2(nrows, &lens, nrows, &starts));
2499     /* loop over new rows determining lens and starting points */
2500     for (i = 0; i < nrows; i++) {
2501       kstart    = ai[irow[i]];
2502       kend      = kstart + ailen[irow[i]];
2503       starts[i] = kstart;
2504       for (k = kstart; k < kend; k++) {
2505         if (aj[k] >= first) {
2506           starts[i] = k;
2507           break;
2508         }
2509       }
2510       sum = 0;
2511       while (k < kend) {
2512         if (aj[k++] >= first + ncols) break;
2513         sum++;
2514       }
2515       lens[i] = sum;
2516     }
2517     /* create submatrix */
2518     if (scall == MAT_REUSE_MATRIX) {
2519       PetscInt n_cols, n_rows;
2520       PetscCall(MatGetSize(*B, &n_rows, &n_cols));
2521       PetscCheck(n_rows == nrows && n_cols == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reused submatrix wrong size");
2522       PetscCall(MatZeroEntries(*B));
2523       C = *B;
2524     } else {
2525       PetscInt rbs, cbs;
2526       PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
2527       PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE));
2528       PetscCall(ISGetBlockSize(isrow, &rbs));
2529       PetscCall(ISGetBlockSize(iscol, &cbs));
2530       PetscCall(MatSetBlockSizes(C, rbs, cbs));
2531       PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
2532       PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(C, 0, lens));
2533     }
2534     c = (Mat_SeqAIJ *)C->data;
2535 
2536     /* loop over rows inserting into submatrix */
2537     PetscCall(MatSeqAIJGetArrayWrite(C, &a_new)); // Not 'a_new = c->a-new', since that raw usage ignores offload state of C
2538     j_new = c->j;
2539     i_new = c->i;
2540     PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2541     for (i = 0; i < nrows; i++) {
2542       ii    = starts[i];
2543       lensi = lens[i];
2544       for (k = 0; k < lensi; k++) *j_new++ = aj[ii + k] - first;
2545       PetscCall(PetscArraycpy(a_new, aa + starts[i], lensi));
2546       a_new += lensi;
2547       i_new[i + 1] = i_new[i] + lensi;
2548       c->ilen[i]   = lensi;
2549     }
2550     PetscCall(MatSeqAIJRestoreArrayWrite(C, &a_new)); // Set C's offload state properly
2551     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2552     PetscCall(PetscFree2(lens, starts));
2553   } else {
2554     PetscCall(ISGetIndices(iscol, &icol));
2555     PetscCall(PetscCalloc1(oldcols, &smap));
2556     PetscCall(PetscMalloc1(1 + nrows, &lens));
2557     for (i = 0; i < ncols; i++) {
2558       PetscCheck(icol[i] < oldcols, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Requesting column beyond largest column icol[%" PetscInt_FMT "] %" PetscInt_FMT " >= A->cmap->n %" PetscInt_FMT, i, icol[i], oldcols);
2559       smap[icol[i]] = i + 1;
2560     }
2561 
2562     /* determine lens of each row */
2563     for (i = 0; i < nrows; i++) {
2564       kstart  = ai[irow[i]];
2565       kend    = kstart + a->ilen[irow[i]];
2566       lens[i] = 0;
2567       for (k = kstart; k < kend; k++) {
2568         if (smap[aj[k]]) lens[i]++;
2569       }
2570     }
2571     /* Create and fill new matrix */
2572     if (scall == MAT_REUSE_MATRIX) {
2573       PetscBool equal;
2574 
2575       c = (Mat_SeqAIJ *)((*B)->data);
2576       PetscCheck((*B)->rmap->n == nrows && (*B)->cmap->n == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong size");
2577       PetscCall(PetscArraycmp(c->ilen, lens, (*B)->rmap->n, &equal));
2578       PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong no of nonzeros");
2579       PetscCall(PetscArrayzero(c->ilen, (*B)->rmap->n));
2580       C = *B;
2581     } else {
2582       PetscInt rbs, cbs;
2583       PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
2584       PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE));
2585       PetscCall(ISGetBlockSize(isrow, &rbs));
2586       PetscCall(ISGetBlockSize(iscol, &cbs));
2587       PetscCall(MatSetBlockSizes(C, rbs, cbs));
2588       PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
2589       PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(C, 0, lens));
2590     }
2591     PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2592 
2593     c = (Mat_SeqAIJ *)(C->data);
2594     PetscCall(MatSeqAIJGetArrayWrite(C, &c_a)); // Not 'c->a', since that raw usage ignores offload state of C
2595     for (i = 0; i < nrows; i++) {
2596       row      = irow[i];
2597       kstart   = ai[row];
2598       kend     = kstart + a->ilen[row];
2599       mat_i    = c->i[i];
2600       mat_j    = c->j + mat_i;
2601       mat_a    = c_a + mat_i;
2602       mat_ilen = c->ilen + i;
2603       for (k = kstart; k < kend; k++) {
2604         if ((tcol = smap[a->j[k]])) {
2605           *mat_j++ = tcol - 1;
2606           *mat_a++ = aa[k];
2607           (*mat_ilen)++;
2608         }
2609       }
2610     }
2611     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2612     /* Free work space */
2613     PetscCall(ISRestoreIndices(iscol, &icol));
2614     PetscCall(PetscFree(smap));
2615     PetscCall(PetscFree(lens));
2616     /* sort */
2617     for (i = 0; i < nrows; i++) {
2618       PetscInt ilen;
2619 
2620       mat_i = c->i[i];
2621       mat_j = c->j + mat_i;
2622       mat_a = c_a + mat_i;
2623       ilen  = c->ilen[i];
2624       PetscCall(PetscSortIntWithScalarArray(ilen, mat_j, mat_a));
2625     }
2626     PetscCall(MatSeqAIJRestoreArrayWrite(C, &c_a));
2627   }
2628 #if defined(PETSC_HAVE_DEVICE)
2629   PetscCall(MatBindToCPU(C, A->boundtocpu));
2630 #endif
2631   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2632   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
2633 
2634   PetscCall(ISRestoreIndices(isrow, &irow));
2635   *B = C;
2636   PetscFunctionReturn(PETSC_SUCCESS);
2637 }
2638 
2639 PetscErrorCode MatGetMultiProcBlock_SeqAIJ(Mat mat, MPI_Comm subComm, MatReuse scall, Mat *subMat)
2640 {
2641   Mat B;
2642 
2643   PetscFunctionBegin;
2644   if (scall == MAT_INITIAL_MATRIX) {
2645     PetscCall(MatCreate(subComm, &B));
2646     PetscCall(MatSetSizes(B, mat->rmap->n, mat->cmap->n, mat->rmap->n, mat->cmap->n));
2647     PetscCall(MatSetBlockSizesFromMats(B, mat, mat));
2648     PetscCall(MatSetType(B, MATSEQAIJ));
2649     PetscCall(MatDuplicateNoCreate_SeqAIJ(B, mat, MAT_COPY_VALUES, PETSC_TRUE));
2650     *subMat = B;
2651   } else {
2652     PetscCall(MatCopy_SeqAIJ(mat, *subMat, SAME_NONZERO_PATTERN));
2653   }
2654   PetscFunctionReturn(PETSC_SUCCESS);
2655 }
2656 
2657 PetscErrorCode MatILUFactor_SeqAIJ(Mat inA, IS row, IS col, const MatFactorInfo *info)
2658 {
2659   Mat_SeqAIJ *a = (Mat_SeqAIJ *)inA->data;
2660   Mat         outA;
2661   PetscBool   row_identity, col_identity;
2662 
2663   PetscFunctionBegin;
2664   PetscCheck(info->levels == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only levels=0 supported for in-place ilu");
2665 
2666   PetscCall(ISIdentity(row, &row_identity));
2667   PetscCall(ISIdentity(col, &col_identity));
2668 
2669   outA             = inA;
2670   outA->factortype = MAT_FACTOR_LU;
2671   PetscCall(PetscFree(inA->solvertype));
2672   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &inA->solvertype));
2673 
2674   PetscCall(PetscObjectReference((PetscObject)row));
2675   PetscCall(ISDestroy(&a->row));
2676 
2677   a->row = row;
2678 
2679   PetscCall(PetscObjectReference((PetscObject)col));
2680   PetscCall(ISDestroy(&a->col));
2681 
2682   a->col = col;
2683 
2684   /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
2685   PetscCall(ISDestroy(&a->icol));
2686   PetscCall(ISInvertPermutation(col, PETSC_DECIDE, &a->icol));
2687 
2688   if (!a->solve_work) { /* this matrix may have been factored before */
2689     PetscCall(PetscMalloc1(inA->rmap->n + 1, &a->solve_work));
2690   }
2691 
2692   PetscCall(MatMarkDiagonal_SeqAIJ(inA));
2693   if (row_identity && col_identity) {
2694     PetscCall(MatLUFactorNumeric_SeqAIJ_inplace(outA, inA, info));
2695   } else {
2696     PetscCall(MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(outA, inA, info));
2697   }
2698   PetscFunctionReturn(PETSC_SUCCESS);
2699 }
2700 
2701 PetscErrorCode MatScale_SeqAIJ(Mat inA, PetscScalar alpha)
2702 {
2703   Mat_SeqAIJ  *a = (Mat_SeqAIJ *)inA->data;
2704   PetscScalar *v;
2705   PetscBLASInt one = 1, bnz;
2706 
2707   PetscFunctionBegin;
2708   PetscCall(MatSeqAIJGetArray(inA, &v));
2709   PetscCall(PetscBLASIntCast(a->nz, &bnz));
2710   PetscCallBLAS("BLASscal", BLASscal_(&bnz, &alpha, v, &one));
2711   PetscCall(PetscLogFlops(a->nz));
2712   PetscCall(MatSeqAIJRestoreArray(inA, &v));
2713   PetscCall(MatSeqAIJInvalidateDiagonal(inA));
2714   PetscFunctionReturn(PETSC_SUCCESS);
2715 }
2716 
2717 PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *submatj)
2718 {
2719   PetscInt i;
2720 
2721   PetscFunctionBegin;
2722   if (!submatj->id) { /* delete data that are linked only to submats[id=0] */
2723     PetscCall(PetscFree4(submatj->sbuf1, submatj->ptr, submatj->tmp, submatj->ctr));
2724 
2725     for (i = 0; i < submatj->nrqr; ++i) PetscCall(PetscFree(submatj->sbuf2[i]));
2726     PetscCall(PetscFree3(submatj->sbuf2, submatj->req_size, submatj->req_source1));
2727 
2728     if (submatj->rbuf1) {
2729       PetscCall(PetscFree(submatj->rbuf1[0]));
2730       PetscCall(PetscFree(submatj->rbuf1));
2731     }
2732 
2733     for (i = 0; i < submatj->nrqs; ++i) PetscCall(PetscFree(submatj->rbuf3[i]));
2734     PetscCall(PetscFree3(submatj->req_source2, submatj->rbuf2, submatj->rbuf3));
2735     PetscCall(PetscFree(submatj->pa));
2736   }
2737 
2738 #if defined(PETSC_USE_CTABLE)
2739   PetscCall(PetscHMapIDestroy(&submatj->rmap));
2740   if (submatj->cmap_loc) PetscCall(PetscFree(submatj->cmap_loc));
2741   PetscCall(PetscFree(submatj->rmap_loc));
2742 #else
2743   PetscCall(PetscFree(submatj->rmap));
2744 #endif
2745 
2746   if (!submatj->allcolumns) {
2747 #if defined(PETSC_USE_CTABLE)
2748     PetscCall(PetscHMapIDestroy((PetscHMapI *)&submatj->cmap));
2749 #else
2750     PetscCall(PetscFree(submatj->cmap));
2751 #endif
2752   }
2753   PetscCall(PetscFree(submatj->row2proc));
2754 
2755   PetscCall(PetscFree(submatj));
2756   PetscFunctionReturn(PETSC_SUCCESS);
2757 }
2758 
2759 PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat C)
2760 {
2761   Mat_SeqAIJ  *c       = (Mat_SeqAIJ *)C->data;
2762   Mat_SubSppt *submatj = c->submatis1;
2763 
2764   PetscFunctionBegin;
2765   PetscCall((*submatj->destroy)(C));
2766   PetscCall(MatDestroySubMatrix_Private(submatj));
2767   PetscFunctionReturn(PETSC_SUCCESS);
2768 }
2769 
2770 /* Note this has code duplication with MatDestroySubMatrices_SeqBAIJ() */
2771 PetscErrorCode MatDestroySubMatrices_SeqAIJ(PetscInt n, Mat *mat[])
2772 {
2773   PetscInt     i;
2774   Mat          C;
2775   Mat_SeqAIJ  *c;
2776   Mat_SubSppt *submatj;
2777 
2778   PetscFunctionBegin;
2779   for (i = 0; i < n; i++) {
2780     C       = (*mat)[i];
2781     c       = (Mat_SeqAIJ *)C->data;
2782     submatj = c->submatis1;
2783     if (submatj) {
2784       if (--((PetscObject)C)->refct <= 0) {
2785         PetscCall(PetscFree(C->factorprefix));
2786         PetscCall((*submatj->destroy)(C));
2787         PetscCall(MatDestroySubMatrix_Private(submatj));
2788         PetscCall(PetscFree(C->defaultvectype));
2789         PetscCall(PetscFree(C->defaultrandtype));
2790         PetscCall(PetscLayoutDestroy(&C->rmap));
2791         PetscCall(PetscLayoutDestroy(&C->cmap));
2792         PetscCall(PetscHeaderDestroy(&C));
2793       }
2794     } else {
2795       PetscCall(MatDestroy(&C));
2796     }
2797   }
2798 
2799   /* Destroy Dummy submatrices created for reuse */
2800   PetscCall(MatDestroySubMatrices_Dummy(n, mat));
2801 
2802   PetscCall(PetscFree(*mat));
2803   PetscFunctionReturn(PETSC_SUCCESS);
2804 }
2805 
2806 PetscErrorCode MatCreateSubMatrices_SeqAIJ(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *B[])
2807 {
2808   PetscInt i;
2809 
2810   PetscFunctionBegin;
2811   if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscCalloc1(n + 1, B));
2812 
2813   for (i = 0; i < n; i++) PetscCall(MatCreateSubMatrix_SeqAIJ(A, irow[i], icol[i], PETSC_DECIDE, scall, &(*B)[i]));
2814   PetscFunctionReturn(PETSC_SUCCESS);
2815 }
2816 
2817 PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A, PetscInt is_max, IS is[], PetscInt ov)
2818 {
2819   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
2820   PetscInt        row, i, j, k, l, ll, m, n, *nidx, isz, val;
2821   const PetscInt *idx;
2822   PetscInt        start, end, *ai, *aj, bs = (A->rmap->bs > 0 && A->rmap->bs == A->cmap->bs) ? A->rmap->bs : 1;
2823   PetscBT         table;
2824 
2825   PetscFunctionBegin;
2826   m  = A->rmap->n / bs;
2827   ai = a->i;
2828   aj = a->j;
2829 
2830   PetscCheck(ov >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "illegal negative overlap value used");
2831 
2832   PetscCall(PetscMalloc1(m + 1, &nidx));
2833   PetscCall(PetscBTCreate(m, &table));
2834 
2835   for (i = 0; i < is_max; i++) {
2836     /* Initialize the two local arrays */
2837     isz = 0;
2838     PetscCall(PetscBTMemzero(m, table));
2839 
2840     /* Extract the indices, assume there can be duplicate entries */
2841     PetscCall(ISGetIndices(is[i], &idx));
2842     PetscCall(ISGetLocalSize(is[i], &n));
2843 
2844     if (bs > 1) {
2845       /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
2846       for (j = 0; j < n; ++j) {
2847         if (!PetscBTLookupSet(table, idx[j] / bs)) nidx[isz++] = idx[j] / bs;
2848       }
2849       PetscCall(ISRestoreIndices(is[i], &idx));
2850       PetscCall(ISDestroy(&is[i]));
2851 
2852       k = 0;
2853       for (j = 0; j < ov; j++) { /* for each overlap */
2854         n = isz;
2855         for (; k < n; k++) { /* do only those rows in nidx[k], which are not done yet */
2856           for (ll = 0; ll < bs; ll++) {
2857             row   = bs * nidx[k] + ll;
2858             start = ai[row];
2859             end   = ai[row + 1];
2860             for (l = start; l < end; l++) {
2861               val = aj[l] / bs;
2862               if (!PetscBTLookupSet(table, val)) nidx[isz++] = val;
2863             }
2864           }
2865         }
2866       }
2867       PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, isz, nidx, PETSC_COPY_VALUES, (is + i)));
2868     } else {
2869       /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
2870       for (j = 0; j < n; ++j) {
2871         if (!PetscBTLookupSet(table, idx[j])) nidx[isz++] = idx[j];
2872       }
2873       PetscCall(ISRestoreIndices(is[i], &idx));
2874       PetscCall(ISDestroy(&is[i]));
2875 
2876       k = 0;
2877       for (j = 0; j < ov; j++) { /* for each overlap */
2878         n = isz;
2879         for (; k < n; k++) { /* do only those rows in nidx[k], which are not done yet */
2880           row   = nidx[k];
2881           start = ai[row];
2882           end   = ai[row + 1];
2883           for (l = start; l < end; l++) {
2884             val = aj[l];
2885             if (!PetscBTLookupSet(table, val)) nidx[isz++] = val;
2886           }
2887         }
2888       }
2889       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, PETSC_COPY_VALUES, (is + i)));
2890     }
2891   }
2892   PetscCall(PetscBTDestroy(&table));
2893   PetscCall(PetscFree(nidx));
2894   PetscFunctionReturn(PETSC_SUCCESS);
2895 }
2896 
2897 PetscErrorCode MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
2898 {
2899   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
2900   PetscInt        i, nz = 0, m = A->rmap->n, n = A->cmap->n;
2901   const PetscInt *row, *col;
2902   PetscInt       *cnew, j, *lens;
2903   IS              icolp, irowp;
2904   PetscInt       *cwork = NULL;
2905   PetscScalar    *vwork = NULL;
2906 
2907   PetscFunctionBegin;
2908   PetscCall(ISInvertPermutation(rowp, PETSC_DECIDE, &irowp));
2909   PetscCall(ISGetIndices(irowp, &row));
2910   PetscCall(ISInvertPermutation(colp, PETSC_DECIDE, &icolp));
2911   PetscCall(ISGetIndices(icolp, &col));
2912 
2913   /* determine lengths of permuted rows */
2914   PetscCall(PetscMalloc1(m + 1, &lens));
2915   for (i = 0; i < m; i++) lens[row[i]] = a->i[i + 1] - a->i[i];
2916   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
2917   PetscCall(MatSetSizes(*B, m, n, m, n));
2918   PetscCall(MatSetBlockSizesFromMats(*B, A, A));
2919   PetscCall(MatSetType(*B, ((PetscObject)A)->type_name));
2920   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*B, 0, lens));
2921   PetscCall(PetscFree(lens));
2922 
2923   PetscCall(PetscMalloc1(n, &cnew));
2924   for (i = 0; i < m; i++) {
2925     PetscCall(MatGetRow_SeqAIJ(A, i, &nz, &cwork, &vwork));
2926     for (j = 0; j < nz; j++) cnew[j] = col[cwork[j]];
2927     PetscCall(MatSetValues_SeqAIJ(*B, 1, &row[i], nz, cnew, vwork, INSERT_VALUES));
2928     PetscCall(MatRestoreRow_SeqAIJ(A, i, &nz, &cwork, &vwork));
2929   }
2930   PetscCall(PetscFree(cnew));
2931 
2932   (*B)->assembled = PETSC_FALSE;
2933 
2934 #if defined(PETSC_HAVE_DEVICE)
2935   PetscCall(MatBindToCPU(*B, A->boundtocpu));
2936 #endif
2937   PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY));
2938   PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY));
2939   PetscCall(ISRestoreIndices(irowp, &row));
2940   PetscCall(ISRestoreIndices(icolp, &col));
2941   PetscCall(ISDestroy(&irowp));
2942   PetscCall(ISDestroy(&icolp));
2943   if (rowp == colp) PetscCall(MatPropagateSymmetryOptions(A, *B));
2944   PetscFunctionReturn(PETSC_SUCCESS);
2945 }
2946 
2947 PetscErrorCode MatCopy_SeqAIJ(Mat A, Mat B, MatStructure str)
2948 {
2949   PetscFunctionBegin;
2950   /* If the two matrices have the same copy implementation, use fast copy. */
2951   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
2952     Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2953     Mat_SeqAIJ        *b = (Mat_SeqAIJ *)B->data;
2954     const PetscScalar *aa;
2955 
2956     PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2957     PetscCheck(a->i[A->rmap->n] == b->i[B->rmap->n], PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of nonzeros in two matrices are different %" PetscInt_FMT " != %" PetscInt_FMT, a->i[A->rmap->n], b->i[B->rmap->n]);
2958     PetscCall(PetscArraycpy(b->a, aa, a->i[A->rmap->n]));
2959     PetscCall(PetscObjectStateIncrease((PetscObject)B));
2960     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2961   } else {
2962     PetscCall(MatCopy_Basic(A, B, str));
2963   }
2964   PetscFunctionReturn(PETSC_SUCCESS);
2965 }
2966 
2967 PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat A, PetscScalar *array[])
2968 {
2969   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
2970 
2971   PetscFunctionBegin;
2972   *array = a->a;
2973   PetscFunctionReturn(PETSC_SUCCESS);
2974 }
2975 
2976 PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat A, PetscScalar *array[])
2977 {
2978   PetscFunctionBegin;
2979   *array = NULL;
2980   PetscFunctionReturn(PETSC_SUCCESS);
2981 }
2982 
2983 /*
2984    Computes the number of nonzeros per row needed for preallocation when X and Y
2985    have different nonzero structure.
2986 */
2987 PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt m, const PetscInt *xi, const PetscInt *xj, const PetscInt *yi, const PetscInt *yj, PetscInt *nnz)
2988 {
2989   PetscInt i, j, k, nzx, nzy;
2990 
2991   PetscFunctionBegin;
2992   /* Set the number of nonzeros in the new matrix */
2993   for (i = 0; i < m; i++) {
2994     const PetscInt *xjj = xj + xi[i], *yjj = yj + yi[i];
2995     nzx    = xi[i + 1] - xi[i];
2996     nzy    = yi[i + 1] - yi[i];
2997     nnz[i] = 0;
2998     for (j = 0, k = 0; j < nzx; j++) {                  /* Point in X */
2999       for (; k < nzy && yjj[k] < xjj[j]; k++) nnz[i]++; /* Catch up to X */
3000       if (k < nzy && yjj[k] == xjj[j]) k++;             /* Skip duplicate */
3001       nnz[i]++;
3002     }
3003     for (; k < nzy; k++) nnz[i]++;
3004   }
3005   PetscFunctionReturn(PETSC_SUCCESS);
3006 }
3007 
3008 PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat Y, Mat X, PetscInt *nnz)
3009 {
3010   PetscInt    m = Y->rmap->N;
3011   Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data;
3012   Mat_SeqAIJ *y = (Mat_SeqAIJ *)Y->data;
3013 
3014   PetscFunctionBegin;
3015   /* Set the number of nonzeros in the new matrix */
3016   PetscCall(MatAXPYGetPreallocation_SeqX_private(m, x->i, x->j, y->i, y->j, nnz));
3017   PetscFunctionReturn(PETSC_SUCCESS);
3018 }
3019 
3020 PetscErrorCode MatAXPY_SeqAIJ(Mat Y, PetscScalar a, Mat X, MatStructure str)
3021 {
3022   Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data, *y = (Mat_SeqAIJ *)Y->data;
3023 
3024   PetscFunctionBegin;
3025   if (str == UNKNOWN_NONZERO_PATTERN || (PetscDefined(USE_DEBUG) && str == SAME_NONZERO_PATTERN)) {
3026     PetscBool e = x->nz == y->nz ? PETSC_TRUE : PETSC_FALSE;
3027     if (e) {
3028       PetscCall(PetscArraycmp(x->i, y->i, Y->rmap->n + 1, &e));
3029       if (e) {
3030         PetscCall(PetscArraycmp(x->j, y->j, y->nz, &e));
3031         if (e) str = SAME_NONZERO_PATTERN;
3032       }
3033     }
3034     if (!e) PetscCheck(str != SAME_NONZERO_PATTERN, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatStructure is not SAME_NONZERO_PATTERN");
3035   }
3036   if (str == SAME_NONZERO_PATTERN) {
3037     const PetscScalar *xa;
3038     PetscScalar       *ya, alpha = a;
3039     PetscBLASInt       one = 1, bnz;
3040 
3041     PetscCall(PetscBLASIntCast(x->nz, &bnz));
3042     PetscCall(MatSeqAIJGetArray(Y, &ya));
3043     PetscCall(MatSeqAIJGetArrayRead(X, &xa));
3044     PetscCallBLAS("BLASaxpy", BLASaxpy_(&bnz, &alpha, xa, &one, ya, &one));
3045     PetscCall(MatSeqAIJRestoreArrayRead(X, &xa));
3046     PetscCall(MatSeqAIJRestoreArray(Y, &ya));
3047     PetscCall(PetscLogFlops(2.0 * bnz));
3048     PetscCall(MatSeqAIJInvalidateDiagonal(Y));
3049     PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3050   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
3051     PetscCall(MatAXPY_Basic(Y, a, X, str));
3052   } else {
3053     Mat       B;
3054     PetscInt *nnz;
3055     PetscCall(PetscMalloc1(Y->rmap->N, &nnz));
3056     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &B));
3057     PetscCall(PetscObjectSetName((PetscObject)B, ((PetscObject)Y)->name));
3058     PetscCall(MatSetLayouts(B, Y->rmap, Y->cmap));
3059     PetscCall(MatSetType(B, ((PetscObject)Y)->type_name));
3060     PetscCall(MatAXPYGetPreallocation_SeqAIJ(Y, X, nnz));
3061     PetscCall(MatSeqAIJSetPreallocation(B, 0, nnz));
3062     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
3063     PetscCall(MatHeaderMerge(Y, &B));
3064     PetscCall(MatSeqAIJCheckInode(Y));
3065     PetscCall(PetscFree(nnz));
3066   }
3067   PetscFunctionReturn(PETSC_SUCCESS);
3068 }
3069 
3070 PETSC_INTERN PetscErrorCode MatConjugate_SeqAIJ(Mat mat)
3071 {
3072 #if defined(PETSC_USE_COMPLEX)
3073   Mat_SeqAIJ  *aij = (Mat_SeqAIJ *)mat->data;
3074   PetscInt     i, nz;
3075   PetscScalar *a;
3076 
3077   PetscFunctionBegin;
3078   nz = aij->nz;
3079   PetscCall(MatSeqAIJGetArray(mat, &a));
3080   for (i = 0; i < nz; i++) a[i] = PetscConj(a[i]);
3081   PetscCall(MatSeqAIJRestoreArray(mat, &a));
3082 #else
3083   PetscFunctionBegin;
3084 #endif
3085   PetscFunctionReturn(PETSC_SUCCESS);
3086 }
3087 
3088 PetscErrorCode MatGetRowMaxAbs_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3089 {
3090   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3091   PetscInt         i, j, m = A->rmap->n, *ai, *aj, ncols, n;
3092   PetscReal        atmp;
3093   PetscScalar     *x;
3094   const MatScalar *aa, *av;
3095 
3096   PetscFunctionBegin;
3097   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3098   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3099   aa = av;
3100   ai = a->i;
3101   aj = a->j;
3102 
3103   PetscCall(VecSet(v, 0.0));
3104   PetscCall(VecGetArrayWrite(v, &x));
3105   PetscCall(VecGetLocalSize(v, &n));
3106   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3107   for (i = 0; i < m; i++) {
3108     ncols = ai[1] - ai[0];
3109     ai++;
3110     for (j = 0; j < ncols; j++) {
3111       atmp = PetscAbsScalar(*aa);
3112       if (PetscAbsScalar(x[i]) < atmp) {
3113         x[i] = atmp;
3114         if (idx) idx[i] = *aj;
3115       }
3116       aa++;
3117       aj++;
3118     }
3119   }
3120   PetscCall(VecRestoreArrayWrite(v, &x));
3121   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3122   PetscFunctionReturn(PETSC_SUCCESS);
3123 }
3124 
3125 PetscErrorCode MatGetRowMax_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3126 {
3127   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3128   PetscInt         i, j, m = A->rmap->n, *ai, *aj, ncols, n;
3129   PetscScalar     *x;
3130   const MatScalar *aa, *av;
3131 
3132   PetscFunctionBegin;
3133   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3134   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3135   aa = av;
3136   ai = a->i;
3137   aj = a->j;
3138 
3139   PetscCall(VecSet(v, 0.0));
3140   PetscCall(VecGetArrayWrite(v, &x));
3141   PetscCall(VecGetLocalSize(v, &n));
3142   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3143   for (i = 0; i < m; i++) {
3144     ncols = ai[1] - ai[0];
3145     ai++;
3146     if (ncols == A->cmap->n) { /* row is dense */
3147       x[i] = *aa;
3148       if (idx) idx[i] = 0;
3149     } else { /* row is sparse so already KNOW maximum is 0.0 or higher */
3150       x[i] = 0.0;
3151       if (idx) {
3152         for (j = 0; j < ncols; j++) { /* find first implicit 0.0 in the row */
3153           if (aj[j] > j) {
3154             idx[i] = j;
3155             break;
3156           }
3157         }
3158         /* in case first implicit 0.0 in the row occurs at ncols-th column */
3159         if (j == ncols && j < A->cmap->n) idx[i] = j;
3160       }
3161     }
3162     for (j = 0; j < ncols; j++) {
3163       if (PetscRealPart(x[i]) < PetscRealPart(*aa)) {
3164         x[i] = *aa;
3165         if (idx) idx[i] = *aj;
3166       }
3167       aa++;
3168       aj++;
3169     }
3170   }
3171   PetscCall(VecRestoreArrayWrite(v, &x));
3172   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3173   PetscFunctionReturn(PETSC_SUCCESS);
3174 }
3175 
3176 PetscErrorCode MatGetRowMinAbs_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3177 {
3178   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3179   PetscInt         i, j, m = A->rmap->n, *ai, *aj, ncols, n;
3180   PetscScalar     *x;
3181   const MatScalar *aa, *av;
3182 
3183   PetscFunctionBegin;
3184   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3185   aa = av;
3186   ai = a->i;
3187   aj = a->j;
3188 
3189   PetscCall(VecSet(v, 0.0));
3190   PetscCall(VecGetArrayWrite(v, &x));
3191   PetscCall(VecGetLocalSize(v, &n));
3192   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector, %" PetscInt_FMT " vs. %" PetscInt_FMT " rows", m, n);
3193   for (i = 0; i < m; i++) {
3194     ncols = ai[1] - ai[0];
3195     ai++;
3196     if (ncols == A->cmap->n) { /* row is dense */
3197       x[i] = *aa;
3198       if (idx) idx[i] = 0;
3199     } else { /* row is sparse so already KNOW minimum is 0.0 or higher */
3200       x[i] = 0.0;
3201       if (idx) { /* find first implicit 0.0 in the row */
3202         for (j = 0; j < ncols; j++) {
3203           if (aj[j] > j) {
3204             idx[i] = j;
3205             break;
3206           }
3207         }
3208         /* in case first implicit 0.0 in the row occurs at ncols-th column */
3209         if (j == ncols && j < A->cmap->n) idx[i] = j;
3210       }
3211     }
3212     for (j = 0; j < ncols; j++) {
3213       if (PetscAbsScalar(x[i]) > PetscAbsScalar(*aa)) {
3214         x[i] = *aa;
3215         if (idx) idx[i] = *aj;
3216       }
3217       aa++;
3218       aj++;
3219     }
3220   }
3221   PetscCall(VecRestoreArrayWrite(v, &x));
3222   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3223   PetscFunctionReturn(PETSC_SUCCESS);
3224 }
3225 
3226 PetscErrorCode MatGetRowMin_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3227 {
3228   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3229   PetscInt         i, j, m = A->rmap->n, ncols, n;
3230   const PetscInt  *ai, *aj;
3231   PetscScalar     *x;
3232   const MatScalar *aa, *av;
3233 
3234   PetscFunctionBegin;
3235   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3236   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3237   aa = av;
3238   ai = a->i;
3239   aj = a->j;
3240 
3241   PetscCall(VecSet(v, 0.0));
3242   PetscCall(VecGetArrayWrite(v, &x));
3243   PetscCall(VecGetLocalSize(v, &n));
3244   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3245   for (i = 0; i < m; i++) {
3246     ncols = ai[1] - ai[0];
3247     ai++;
3248     if (ncols == A->cmap->n) { /* row is dense */
3249       x[i] = *aa;
3250       if (idx) idx[i] = 0;
3251     } else { /* row is sparse so already KNOW minimum is 0.0 or lower */
3252       x[i] = 0.0;
3253       if (idx) { /* find first implicit 0.0 in the row */
3254         for (j = 0; j < ncols; j++) {
3255           if (aj[j] > j) {
3256             idx[i] = j;
3257             break;
3258           }
3259         }
3260         /* in case first implicit 0.0 in the row occurs at ncols-th column */
3261         if (j == ncols && j < A->cmap->n) idx[i] = j;
3262       }
3263     }
3264     for (j = 0; j < ncols; j++) {
3265       if (PetscRealPart(x[i]) > PetscRealPart(*aa)) {
3266         x[i] = *aa;
3267         if (idx) idx[i] = *aj;
3268       }
3269       aa++;
3270       aj++;
3271     }
3272   }
3273   PetscCall(VecRestoreArrayWrite(v, &x));
3274   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3275   PetscFunctionReturn(PETSC_SUCCESS);
3276 }
3277 
3278 PetscErrorCode MatInvertBlockDiagonal_SeqAIJ(Mat A, const PetscScalar **values)
3279 {
3280   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
3281   PetscInt        i, bs = PetscAbs(A->rmap->bs), mbs = A->rmap->n / bs, ipvt[5], bs2 = bs * bs, *v_pivots, ij[7], *IJ, j;
3282   MatScalar      *diag, work[25], *v_work;
3283   const PetscReal shift = 0.0;
3284   PetscBool       allowzeropivot, zeropivotdetected = PETSC_FALSE;
3285 
3286   PetscFunctionBegin;
3287   allowzeropivot = PetscNot(A->erroriffailure);
3288   if (a->ibdiagvalid) {
3289     if (values) *values = a->ibdiag;
3290     PetscFunctionReturn(PETSC_SUCCESS);
3291   }
3292   PetscCall(MatMarkDiagonal_SeqAIJ(A));
3293   if (!a->ibdiag) { PetscCall(PetscMalloc1(bs2 * mbs, &a->ibdiag)); }
3294   diag = a->ibdiag;
3295   if (values) *values = a->ibdiag;
3296   /* factor and invert each block */
3297   switch (bs) {
3298   case 1:
3299     for (i = 0; i < mbs; i++) {
3300       PetscCall(MatGetValues(A, 1, &i, 1, &i, diag + i));
3301       if (PetscAbsScalar(diag[i] + shift) < PETSC_MACHINE_EPSILON) {
3302         if (allowzeropivot) {
3303           A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3304           A->factorerror_zeropivot_value = PetscAbsScalar(diag[i]);
3305           A->factorerror_zeropivot_row   = i;
3306           PetscCall(PetscInfo(A, "Zero pivot, row %" PetscInt_FMT " pivot %g tolerance %g\n", i, (double)PetscAbsScalar(diag[i]), (double)PETSC_MACHINE_EPSILON));
3307         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT " pivot %g tolerance %g", i, (double)PetscAbsScalar(diag[i]), (double)PETSC_MACHINE_EPSILON);
3308       }
3309       diag[i] = (PetscScalar)1.0 / (diag[i] + shift);
3310     }
3311     break;
3312   case 2:
3313     for (i = 0; i < mbs; i++) {
3314       ij[0] = 2 * i;
3315       ij[1] = 2 * i + 1;
3316       PetscCall(MatGetValues(A, 2, ij, 2, ij, diag));
3317       PetscCall(PetscKernel_A_gets_inverse_A_2(diag, shift, allowzeropivot, &zeropivotdetected));
3318       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3319       PetscCall(PetscKernel_A_gets_transpose_A_2(diag));
3320       diag += 4;
3321     }
3322     break;
3323   case 3:
3324     for (i = 0; i < mbs; i++) {
3325       ij[0] = 3 * i;
3326       ij[1] = 3 * i + 1;
3327       ij[2] = 3 * i + 2;
3328       PetscCall(MatGetValues(A, 3, ij, 3, ij, diag));
3329       PetscCall(PetscKernel_A_gets_inverse_A_3(diag, shift, allowzeropivot, &zeropivotdetected));
3330       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3331       PetscCall(PetscKernel_A_gets_transpose_A_3(diag));
3332       diag += 9;
3333     }
3334     break;
3335   case 4:
3336     for (i = 0; i < mbs; i++) {
3337       ij[0] = 4 * i;
3338       ij[1] = 4 * i + 1;
3339       ij[2] = 4 * i + 2;
3340       ij[3] = 4 * i + 3;
3341       PetscCall(MatGetValues(A, 4, ij, 4, ij, diag));
3342       PetscCall(PetscKernel_A_gets_inverse_A_4(diag, shift, allowzeropivot, &zeropivotdetected));
3343       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3344       PetscCall(PetscKernel_A_gets_transpose_A_4(diag));
3345       diag += 16;
3346     }
3347     break;
3348   case 5:
3349     for (i = 0; i < mbs; i++) {
3350       ij[0] = 5 * i;
3351       ij[1] = 5 * i + 1;
3352       ij[2] = 5 * i + 2;
3353       ij[3] = 5 * i + 3;
3354       ij[4] = 5 * i + 4;
3355       PetscCall(MatGetValues(A, 5, ij, 5, ij, diag));
3356       PetscCall(PetscKernel_A_gets_inverse_A_5(diag, ipvt, work, shift, allowzeropivot, &zeropivotdetected));
3357       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3358       PetscCall(PetscKernel_A_gets_transpose_A_5(diag));
3359       diag += 25;
3360     }
3361     break;
3362   case 6:
3363     for (i = 0; i < mbs; i++) {
3364       ij[0] = 6 * i;
3365       ij[1] = 6 * i + 1;
3366       ij[2] = 6 * i + 2;
3367       ij[3] = 6 * i + 3;
3368       ij[4] = 6 * i + 4;
3369       ij[5] = 6 * i + 5;
3370       PetscCall(MatGetValues(A, 6, ij, 6, ij, diag));
3371       PetscCall(PetscKernel_A_gets_inverse_A_6(diag, shift, allowzeropivot, &zeropivotdetected));
3372       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3373       PetscCall(PetscKernel_A_gets_transpose_A_6(diag));
3374       diag += 36;
3375     }
3376     break;
3377   case 7:
3378     for (i = 0; i < mbs; i++) {
3379       ij[0] = 7 * i;
3380       ij[1] = 7 * i + 1;
3381       ij[2] = 7 * i + 2;
3382       ij[3] = 7 * i + 3;
3383       ij[4] = 7 * i + 4;
3384       ij[5] = 7 * i + 5;
3385       ij[6] = 7 * i + 6;
3386       PetscCall(MatGetValues(A, 7, ij, 7, ij, diag));
3387       PetscCall(PetscKernel_A_gets_inverse_A_7(diag, shift, allowzeropivot, &zeropivotdetected));
3388       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3389       PetscCall(PetscKernel_A_gets_transpose_A_7(diag));
3390       diag += 49;
3391     }
3392     break;
3393   default:
3394     PetscCall(PetscMalloc3(bs, &v_work, bs, &v_pivots, bs, &IJ));
3395     for (i = 0; i < mbs; i++) {
3396       for (j = 0; j < bs; j++) IJ[j] = bs * i + j;
3397       PetscCall(MatGetValues(A, bs, IJ, bs, IJ, diag));
3398       PetscCall(PetscKernel_A_gets_inverse_A(bs, diag, v_pivots, v_work, allowzeropivot, &zeropivotdetected));
3399       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3400       PetscCall(PetscKernel_A_gets_transpose_A_N(diag, bs));
3401       diag += bs2;
3402     }
3403     PetscCall(PetscFree3(v_work, v_pivots, IJ));
3404   }
3405   a->ibdiagvalid = PETSC_TRUE;
3406   PetscFunctionReturn(PETSC_SUCCESS);
3407 }
3408 
3409 static PetscErrorCode MatSetRandom_SeqAIJ(Mat x, PetscRandom rctx)
3410 {
3411   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)x->data;
3412   PetscScalar a, *aa;
3413   PetscInt    m, n, i, j, col;
3414 
3415   PetscFunctionBegin;
3416   if (!x->assembled) {
3417     PetscCall(MatGetSize(x, &m, &n));
3418     for (i = 0; i < m; i++) {
3419       for (j = 0; j < aij->imax[i]; j++) {
3420         PetscCall(PetscRandomGetValue(rctx, &a));
3421         col = (PetscInt)(n * PetscRealPart(a));
3422         PetscCall(MatSetValues(x, 1, &i, 1, &col, &a, ADD_VALUES));
3423       }
3424     }
3425   } else {
3426     PetscCall(MatSeqAIJGetArrayWrite(x, &aa));
3427     for (i = 0; i < aij->nz; i++) PetscCall(PetscRandomGetValue(rctx, aa + i));
3428     PetscCall(MatSeqAIJRestoreArrayWrite(x, &aa));
3429   }
3430   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
3431   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
3432   PetscFunctionReturn(PETSC_SUCCESS);
3433 }
3434 
3435 /* Like MatSetRandom_SeqAIJ, but do not set values on columns in range of [low, high) */
3436 PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat x, PetscInt low, PetscInt high, PetscRandom rctx)
3437 {
3438   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)x->data;
3439   PetscScalar a;
3440   PetscInt    m, n, i, j, col, nskip;
3441 
3442   PetscFunctionBegin;
3443   nskip = high - low;
3444   PetscCall(MatGetSize(x, &m, &n));
3445   n -= nskip; /* shrink number of columns where nonzeros can be set */
3446   for (i = 0; i < m; i++) {
3447     for (j = 0; j < aij->imax[i]; j++) {
3448       PetscCall(PetscRandomGetValue(rctx, &a));
3449       col = (PetscInt)(n * PetscRealPart(a));
3450       if (col >= low) col += nskip; /* shift col rightward to skip the hole */
3451       PetscCall(MatSetValues(x, 1, &i, 1, &col, &a, ADD_VALUES));
3452     }
3453   }
3454   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
3455   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
3456   PetscFunctionReturn(PETSC_SUCCESS);
3457 }
3458 
3459 static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
3460                                        MatGetRow_SeqAIJ,
3461                                        MatRestoreRow_SeqAIJ,
3462                                        MatMult_SeqAIJ,
3463                                        /*  4*/ MatMultAdd_SeqAIJ,
3464                                        MatMultTranspose_SeqAIJ,
3465                                        MatMultTransposeAdd_SeqAIJ,
3466                                        NULL,
3467                                        NULL,
3468                                        NULL,
3469                                        /* 10*/ NULL,
3470                                        MatLUFactor_SeqAIJ,
3471                                        NULL,
3472                                        MatSOR_SeqAIJ,
3473                                        MatTranspose_SeqAIJ,
3474                                        /*1 5*/ MatGetInfo_SeqAIJ,
3475                                        MatEqual_SeqAIJ,
3476                                        MatGetDiagonal_SeqAIJ,
3477                                        MatDiagonalScale_SeqAIJ,
3478                                        MatNorm_SeqAIJ,
3479                                        /* 20*/ NULL,
3480                                        MatAssemblyEnd_SeqAIJ,
3481                                        MatSetOption_SeqAIJ,
3482                                        MatZeroEntries_SeqAIJ,
3483                                        /* 24*/ MatZeroRows_SeqAIJ,
3484                                        NULL,
3485                                        NULL,
3486                                        NULL,
3487                                        NULL,
3488                                        /* 29*/ MatSetUp_Seq_Hash,
3489                                        NULL,
3490                                        NULL,
3491                                        NULL,
3492                                        NULL,
3493                                        /* 34*/ MatDuplicate_SeqAIJ,
3494                                        NULL,
3495                                        NULL,
3496                                        MatILUFactor_SeqAIJ,
3497                                        NULL,
3498                                        /* 39*/ MatAXPY_SeqAIJ,
3499                                        MatCreateSubMatrices_SeqAIJ,
3500                                        MatIncreaseOverlap_SeqAIJ,
3501                                        MatGetValues_SeqAIJ,
3502                                        MatCopy_SeqAIJ,
3503                                        /* 44*/ MatGetRowMax_SeqAIJ,
3504                                        MatScale_SeqAIJ,
3505                                        MatShift_SeqAIJ,
3506                                        MatDiagonalSet_SeqAIJ,
3507                                        MatZeroRowsColumns_SeqAIJ,
3508                                        /* 49*/ MatSetRandom_SeqAIJ,
3509                                        MatGetRowIJ_SeqAIJ,
3510                                        MatRestoreRowIJ_SeqAIJ,
3511                                        MatGetColumnIJ_SeqAIJ,
3512                                        MatRestoreColumnIJ_SeqAIJ,
3513                                        /* 54*/ MatFDColoringCreate_SeqXAIJ,
3514                                        NULL,
3515                                        NULL,
3516                                        MatPermute_SeqAIJ,
3517                                        NULL,
3518                                        /* 59*/ NULL,
3519                                        MatDestroy_SeqAIJ,
3520                                        MatView_SeqAIJ,
3521                                        NULL,
3522                                        NULL,
3523                                        /* 64*/ NULL,
3524                                        MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ,
3525                                        NULL,
3526                                        NULL,
3527                                        NULL,
3528                                        /* 69*/ MatGetRowMaxAbs_SeqAIJ,
3529                                        MatGetRowMinAbs_SeqAIJ,
3530                                        NULL,
3531                                        NULL,
3532                                        NULL,
3533                                        /* 74*/ NULL,
3534                                        MatFDColoringApply_AIJ,
3535                                        NULL,
3536                                        NULL,
3537                                        NULL,
3538                                        /* 79*/ MatFindZeroDiagonals_SeqAIJ,
3539                                        NULL,
3540                                        NULL,
3541                                        NULL,
3542                                        MatLoad_SeqAIJ,
3543                                        /* 84*/ MatIsSymmetric_SeqAIJ,
3544                                        MatIsHermitian_SeqAIJ,
3545                                        NULL,
3546                                        NULL,
3547                                        NULL,
3548                                        /* 89*/ NULL,
3549                                        NULL,
3550                                        MatMatMultNumeric_SeqAIJ_SeqAIJ,
3551                                        NULL,
3552                                        NULL,
3553                                        /* 94*/ MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy,
3554                                        NULL,
3555                                        NULL,
3556                                        MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ,
3557                                        NULL,
3558                                        /* 99*/ MatProductSetFromOptions_SeqAIJ,
3559                                        NULL,
3560                                        NULL,
3561                                        MatConjugate_SeqAIJ,
3562                                        NULL,
3563                                        /*104*/ MatSetValuesRow_SeqAIJ,
3564                                        MatRealPart_SeqAIJ,
3565                                        MatImaginaryPart_SeqAIJ,
3566                                        NULL,
3567                                        NULL,
3568                                        /*109*/ MatMatSolve_SeqAIJ,
3569                                        NULL,
3570                                        MatGetRowMin_SeqAIJ,
3571                                        NULL,
3572                                        MatMissingDiagonal_SeqAIJ,
3573                                        /*114*/ NULL,
3574                                        NULL,
3575                                        NULL,
3576                                        NULL,
3577                                        NULL,
3578                                        /*119*/ NULL,
3579                                        NULL,
3580                                        NULL,
3581                                        NULL,
3582                                        MatGetMultiProcBlock_SeqAIJ,
3583                                        /*124*/ MatFindNonzeroRows_SeqAIJ,
3584                                        MatGetColumnReductions_SeqAIJ,
3585                                        MatInvertBlockDiagonal_SeqAIJ,
3586                                        MatInvertVariableBlockDiagonal_SeqAIJ,
3587                                        NULL,
3588                                        /*129*/ NULL,
3589                                        NULL,
3590                                        NULL,
3591                                        MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ,
3592                                        MatTransposeColoringCreate_SeqAIJ,
3593                                        /*134*/ MatTransColoringApplySpToDen_SeqAIJ,
3594                                        MatTransColoringApplyDenToSp_SeqAIJ,
3595                                        NULL,
3596                                        NULL,
3597                                        MatRARtNumeric_SeqAIJ_SeqAIJ,
3598                                        /*139*/ NULL,
3599                                        NULL,
3600                                        NULL,
3601                                        MatFDColoringSetUp_SeqXAIJ,
3602                                        MatFindOffBlockDiagonalEntries_SeqAIJ,
3603                                        MatCreateMPIMatConcatenateSeqMat_SeqAIJ,
3604                                        /*145*/ MatDestroySubMatrices_SeqAIJ,
3605                                        NULL,
3606                                        NULL,
3607                                        MatCreateGraph_Simple_AIJ,
3608                                        NULL,
3609                                        /*150*/ MatTransposeSymbolic_SeqAIJ,
3610                                        MatEliminateZeros_SeqAIJ};
3611 
3612 PetscErrorCode MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat, PetscInt *indices)
3613 {
3614   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
3615   PetscInt    i, nz, n;
3616 
3617   PetscFunctionBegin;
3618   nz = aij->maxnz;
3619   n  = mat->rmap->n;
3620   for (i = 0; i < nz; i++) aij->j[i] = indices[i];
3621   aij->nz = nz;
3622   for (i = 0; i < n; i++) aij->ilen[i] = aij->imax[i];
3623   PetscFunctionReturn(PETSC_SUCCESS);
3624 }
3625 
3626 /*
3627  * Given a sparse matrix with global column indices, compact it by using a local column space.
3628  * The result matrix helps saving memory in other algorithms, such as MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable()
3629  */
3630 PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat mat, ISLocalToGlobalMapping *mapping)
3631 {
3632   Mat_SeqAIJ   *aij = (Mat_SeqAIJ *)mat->data;
3633   PetscHMapI    gid1_lid1;
3634   PetscHashIter tpos;
3635   PetscInt      gid, lid, i, ec, nz = aij->nz;
3636   PetscInt     *garray, *jj = aij->j;
3637 
3638   PetscFunctionBegin;
3639   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
3640   PetscValidPointer(mapping, 2);
3641   /* use a table */
3642   PetscCall(PetscHMapICreateWithSize(mat->rmap->n, &gid1_lid1));
3643   ec = 0;
3644   for (i = 0; i < nz; i++) {
3645     PetscInt data, gid1 = jj[i] + 1;
3646     PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &data));
3647     if (!data) {
3648       /* one based table */
3649       PetscCall(PetscHMapISet(gid1_lid1, gid1, ++ec));
3650     }
3651   }
3652   /* form array of columns we need */
3653   PetscCall(PetscMalloc1(ec, &garray));
3654   PetscHashIterBegin(gid1_lid1, tpos);
3655   while (!PetscHashIterAtEnd(gid1_lid1, tpos)) {
3656     PetscHashIterGetKey(gid1_lid1, tpos, gid);
3657     PetscHashIterGetVal(gid1_lid1, tpos, lid);
3658     PetscHashIterNext(gid1_lid1, tpos);
3659     gid--;
3660     lid--;
3661     garray[lid] = gid;
3662   }
3663   PetscCall(PetscSortInt(ec, garray)); /* sort, and rebuild */
3664   PetscCall(PetscHMapIClear(gid1_lid1));
3665   for (i = 0; i < ec; i++) PetscCall(PetscHMapISet(gid1_lid1, garray[i] + 1, i + 1));
3666   /* compact out the extra columns in B */
3667   for (i = 0; i < nz; i++) {
3668     PetscInt gid1 = jj[i] + 1;
3669     PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &lid));
3670     lid--;
3671     jj[i] = lid;
3672   }
3673   PetscCall(PetscLayoutDestroy(&mat->cmap));
3674   PetscCall(PetscHMapIDestroy(&gid1_lid1));
3675   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)mat), ec, ec, 1, &mat->cmap));
3676   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, mat->cmap->bs, mat->cmap->n, garray, PETSC_OWN_POINTER, mapping));
3677   PetscCall(ISLocalToGlobalMappingSetType(*mapping, ISLOCALTOGLOBALMAPPINGHASH));
3678   PetscFunctionReturn(PETSC_SUCCESS);
3679 }
3680 
3681 /*@
3682     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
3683        in the matrix.
3684 
3685   Input Parameters:
3686 +  mat - the `MATSEQAIJ` matrix
3687 -  indices - the column indices
3688 
3689   Level: advanced
3690 
3691   Notes:
3692     This can be called if you have precomputed the nonzero structure of the
3693   matrix and want to provide it to the matrix object to improve the performance
3694   of the `MatSetValues()` operation.
3695 
3696     You MUST have set the correct numbers of nonzeros per row in the call to
3697   `MatCreateSeqAIJ()`, and the columns indices MUST be sorted.
3698 
3699     MUST be called before any calls to `MatSetValues()`
3700 
3701     The indices should start with zero, not one.
3702 
3703 .seealso: [](ch_matrices), `Mat`, `MATSEQAIJ`
3704 @*/
3705 PetscErrorCode MatSeqAIJSetColumnIndices(Mat mat, PetscInt *indices)
3706 {
3707   PetscFunctionBegin;
3708   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
3709   PetscValidIntPointer(indices, 2);
3710   PetscUseMethod(mat, "MatSeqAIJSetColumnIndices_C", (Mat, PetscInt *), (mat, indices));
3711   PetscFunctionReturn(PETSC_SUCCESS);
3712 }
3713 
3714 PetscErrorCode MatStoreValues_SeqAIJ(Mat mat)
3715 {
3716   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
3717   size_t      nz  = aij->i[mat->rmap->n];
3718 
3719   PetscFunctionBegin;
3720   PetscCheck(aij->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
3721 
3722   /* allocate space for values if not already there */
3723   if (!aij->saved_values) { PetscCall(PetscMalloc1(nz + 1, &aij->saved_values)); }
3724 
3725   /* copy values over */
3726   PetscCall(PetscArraycpy(aij->saved_values, aij->a, nz));
3727   PetscFunctionReturn(PETSC_SUCCESS);
3728 }
3729 
3730 /*@
3731     MatStoreValues - Stashes a copy of the matrix values; this allows reusing of the linear part of a Jacobian, while recomputing only the
3732     nonlinear portion.
3733 
3734    Logically Collect
3735 
3736   Input Parameter:
3737 .  mat - the matrix (currently only `MATAIJ` matrices support this option)
3738 
3739   Level: advanced
3740 
3741   Usage:
3742 .vb
3743     Using SNES
3744     Create Jacobian matrix
3745     Set linear terms into matrix
3746     Apply boundary conditions to matrix, at this time matrix must have
3747       final nonzero structure (i.e. setting the nonlinear terms and applying
3748       boundary conditions again will not change the nonzero structure
3749     MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
3750     MatStoreValues(mat);
3751     Call SNESSetJacobian() with matrix
3752     In your Jacobian routine
3753       MatRetrieveValues(mat);
3754       Set nonlinear terms in matrix
3755 
3756     Without `SNESSolve()`, i.e. when you handle nonlinear solve yourself:
3757     // build linear portion of Jacobian
3758     MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
3759     MatStoreValues(mat);
3760     loop over nonlinear iterations
3761        MatRetrieveValues(mat);
3762        // call MatSetValues(mat,...) to set nonliner portion of Jacobian
3763        // call MatAssemblyBegin/End() on matrix
3764        Solve linear system with Jacobian
3765     endloop
3766 .ve
3767 
3768   Notes:
3769     Matrix must already be assembled before calling this routine
3770     Must set the matrix option `MatSetOption`(mat,`MAT_NEW_NONZERO_LOCATIONS`,`PETSC_FALSE`); before
3771     calling this routine.
3772 
3773     When this is called multiple times it overwrites the previous set of stored values
3774     and does not allocated additional space.
3775 
3776 .seealso: [](ch_matrices), `Mat`, `Mat`, `MatRetrieveValues()`
3777 @*/
3778 PetscErrorCode MatStoreValues(Mat mat)
3779 {
3780   PetscFunctionBegin;
3781   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
3782   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3783   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3784   PetscUseMethod(mat, "MatStoreValues_C", (Mat), (mat));
3785   PetscFunctionReturn(PETSC_SUCCESS);
3786 }
3787 
3788 PetscErrorCode MatRetrieveValues_SeqAIJ(Mat mat)
3789 {
3790   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
3791   PetscInt    nz  = aij->i[mat->rmap->n];
3792 
3793   PetscFunctionBegin;
3794   PetscCheck(aij->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
3795   PetscCheck(aij->saved_values, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatStoreValues(A);first");
3796   /* copy values over */
3797   PetscCall(PetscArraycpy(aij->a, aij->saved_values, nz));
3798   PetscFunctionReturn(PETSC_SUCCESS);
3799 }
3800 
3801 /*@
3802     MatRetrieveValues - Retrieves the copy of the matrix values that was stored with `MatStoreValues()`
3803 
3804    Logically Collect
3805 
3806   Input Parameter:
3807 .  mat - the matrix (currently only `MATAIJ` matrices support this option)
3808 
3809   Level: advanced
3810 
3811 .seealso: [](ch_matrices), `Mat`, `MatStoreValues()`
3812 @*/
3813 PetscErrorCode MatRetrieveValues(Mat mat)
3814 {
3815   PetscFunctionBegin;
3816   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
3817   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3818   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3819   PetscUseMethod(mat, "MatRetrieveValues_C", (Mat), (mat));
3820   PetscFunctionReturn(PETSC_SUCCESS);
3821 }
3822 
3823 /*@C
3824    MatCreateSeqAIJ - Creates a sparse matrix in `MATSEQAIJ` (compressed row) format
3825    (the default parallel PETSc format).  For good matrix assembly performance
3826    the user should preallocate the matrix storage by setting the parameter `nz`
3827    (or the array `nnz`).
3828 
3829    Collective
3830 
3831    Input Parameters:
3832 +  comm - MPI communicator, set to `PETSC_COMM_SELF`
3833 .  m - number of rows
3834 .  n - number of columns
3835 .  nz - number of nonzeros per row (same for all rows)
3836 -  nnz - array containing the number of nonzeros in the various rows
3837          (possibly different for each row) or NULL
3838 
3839    Output Parameter:
3840 .  A - the matrix
3841 
3842    Options Database Keys:
3843 +  -mat_no_inode  - Do not use inodes
3844 -  -mat_inode_limit <limit> - Sets inode limit (max limit=5)
3845 
3846    Level: intermediate
3847 
3848    Notes:
3849    If `nnz` is given then `nz` is ignored
3850 
3851    The `MATSEQAIJ` format, also called
3852    compressed row storage, is fully compatible with standard Fortran
3853    storage.  That is, the stored row and column indices can begin at
3854    either one (as in Fortran) or zero.
3855 
3856    Specify the preallocated storage with either `nz` or `nnz` (not both).
3857    Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
3858    allocation.
3859 
3860    By default, this format uses inodes (identical nodes) when possible, to
3861    improve numerical efficiency of matrix-vector products and solves. We
3862    search for consecutive rows with the same nonzero structure, thereby
3863    reusing matrix information to achieve increased efficiency.
3864 
3865 .seealso: [](ch_matrices), `Mat`, [Sparse Matrix Creation](sec_matsparse), `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()`
3866 @*/
3867 PetscErrorCode MatCreateSeqAIJ(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A)
3868 {
3869   PetscFunctionBegin;
3870   PetscCall(MatCreate(comm, A));
3871   PetscCall(MatSetSizes(*A, m, n, m, n));
3872   PetscCall(MatSetType(*A, MATSEQAIJ));
3873   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*A, nz, nnz));
3874   PetscFunctionReturn(PETSC_SUCCESS);
3875 }
3876 
3877 /*@C
3878    MatSeqAIJSetPreallocation - For good matrix assembly performance
3879    the user should preallocate the matrix storage by setting the parameter nz
3880    (or the array nnz).  By setting these parameters accurately, performance
3881    during matrix assembly can be increased by more than a factor of 50.
3882 
3883    Collective
3884 
3885    Input Parameters:
3886 +  B - The matrix
3887 .  nz - number of nonzeros per row (same for all rows)
3888 -  nnz - array containing the number of nonzeros in the various rows
3889          (possibly different for each row) or NULL
3890 
3891    Options Database Keys:
3892 +  -mat_no_inode  - Do not use inodes
3893 -  -mat_inode_limit <limit> - Sets inode limit (max limit=5)
3894 
3895    Level: intermediate
3896 
3897    Notes:
3898      If `nnz` is given then `nz` is ignored
3899 
3900     The `MATSEQAIJ` format also called
3901    compressed row storage, is fully compatible with standard Fortran
3902    storage.  That is, the stored row and column indices can begin at
3903    either one (as in Fortran) or zero.  See the users' manual for details.
3904 
3905    Specify the preallocated storage with either `nz` or `nnz` (not both).
3906    Set nz = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
3907    allocation.
3908 
3909    You can call `MatGetInfo()` to get information on how effective the preallocation was;
3910    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
3911    You can also run with the option -info and look for messages with the string
3912    malloc in them to see if additional memory allocation was needed.
3913 
3914    Developer Notes:
3915    Use nz of `MAT_SKIP_ALLOCATION` to not allocate any space for the matrix
3916    entries or columns indices
3917 
3918    By default, this format uses inodes (identical nodes) when possible, to
3919    improve numerical efficiency of matrix-vector products and solves. We
3920    search for consecutive rows with the same nonzero structure, thereby
3921    reusing matrix information to achieve increased efficiency.
3922 
3923 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()`, `MatGetInfo()`,
3924           `MatSeqAIJSetTotalPreallocation()`
3925 @*/
3926 PetscErrorCode MatSeqAIJSetPreallocation(Mat B, PetscInt nz, const PetscInt nnz[])
3927 {
3928   PetscFunctionBegin;
3929   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
3930   PetscValidType(B, 1);
3931   PetscTryMethod(B, "MatSeqAIJSetPreallocation_C", (Mat, PetscInt, const PetscInt[]), (B, nz, nnz));
3932   PetscFunctionReturn(PETSC_SUCCESS);
3933 }
3934 
3935 PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat B, PetscInt nz, const PetscInt *nnz)
3936 {
3937   Mat_SeqAIJ *b              = (Mat_SeqAIJ *)B->data;
3938   PetscBool   skipallocation = PETSC_FALSE, realalloc = PETSC_FALSE;
3939   PetscInt    i;
3940 
3941   PetscFunctionBegin;
3942   if (B->hash_active) {
3943     B->ops[0] = b->cops;
3944     PetscCall(PetscHMapIJVDestroy(&b->ht));
3945     PetscCall(PetscFree(b->dnz));
3946     B->hash_active = PETSC_FALSE;
3947   }
3948   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
3949   if (nz == MAT_SKIP_ALLOCATION) {
3950     skipallocation = PETSC_TRUE;
3951     nz             = 0;
3952   }
3953   PetscCall(PetscLayoutSetUp(B->rmap));
3954   PetscCall(PetscLayoutSetUp(B->cmap));
3955 
3956   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
3957   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nz cannot be less than 0: value %" PetscInt_FMT, nz);
3958   if (PetscUnlikelyDebug(nnz)) {
3959     for (i = 0; i < B->rmap->n; i++) {
3960       PetscCheck(nnz[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nnz cannot be less than 0: local row %" PetscInt_FMT " value %" PetscInt_FMT, i, nnz[i]);
3961       PetscCheck(nnz[i] <= B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nnz cannot be greater than row length: local row %" PetscInt_FMT " value %" PetscInt_FMT " rowlength %" PetscInt_FMT, i, nnz[i], B->cmap->n);
3962     }
3963   }
3964 
3965   B->preallocated = PETSC_TRUE;
3966   if (!skipallocation) {
3967     if (!b->imax) { PetscCall(PetscMalloc1(B->rmap->n, &b->imax)); }
3968     if (!b->ilen) {
3969       /* b->ilen will count nonzeros in each row so far. */
3970       PetscCall(PetscCalloc1(B->rmap->n, &b->ilen));
3971     } else {
3972       PetscCall(PetscMemzero(b->ilen, B->rmap->n * sizeof(PetscInt)));
3973     }
3974     if (!b->ipre) PetscCall(PetscMalloc1(B->rmap->n, &b->ipre));
3975     if (!nnz) {
3976       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
3977       else if (nz < 0) nz = 1;
3978       nz = PetscMin(nz, B->cmap->n);
3979       for (i = 0; i < B->rmap->n; i++) b->imax[i] = nz;
3980       nz = nz * B->rmap->n;
3981     } else {
3982       PetscInt64 nz64 = 0;
3983       for (i = 0; i < B->rmap->n; i++) {
3984         b->imax[i] = nnz[i];
3985         nz64 += nnz[i];
3986       }
3987       PetscCall(PetscIntCast(nz64, &nz));
3988     }
3989 
3990     /* allocate the matrix space */
3991     /* FIXME: should B's old memory be unlogged? */
3992     PetscCall(MatSeqXAIJFreeAIJ(B, &b->a, &b->j, &b->i));
3993     if (B->structure_only) {
3994       PetscCall(PetscMalloc1(nz, &b->j));
3995       PetscCall(PetscMalloc1(B->rmap->n + 1, &b->i));
3996     } else {
3997       PetscCall(PetscMalloc3(nz, &b->a, nz, &b->j, B->rmap->n + 1, &b->i));
3998     }
3999     b->i[0] = 0;
4000     for (i = 1; i < B->rmap->n + 1; i++) b->i[i] = b->i[i - 1] + b->imax[i - 1];
4001     if (B->structure_only) {
4002       b->singlemalloc = PETSC_FALSE;
4003       b->free_a       = PETSC_FALSE;
4004     } else {
4005       b->singlemalloc = PETSC_TRUE;
4006       b->free_a       = PETSC_TRUE;
4007     }
4008     b->free_ij = PETSC_TRUE;
4009   } else {
4010     b->free_a  = PETSC_FALSE;
4011     b->free_ij = PETSC_FALSE;
4012   }
4013 
4014   if (b->ipre && nnz != b->ipre && b->imax) {
4015     /* reserve user-requested sparsity */
4016     PetscCall(PetscArraycpy(b->ipre, b->imax, B->rmap->n));
4017   }
4018 
4019   b->nz               = 0;
4020   b->maxnz            = nz;
4021   B->info.nz_unneeded = (double)b->maxnz;
4022   if (realalloc) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
4023   B->was_assembled = PETSC_FALSE;
4024   B->assembled     = PETSC_FALSE;
4025   /* We simply deem preallocation has changed nonzero state. Updating the state
4026      will give clients (like AIJKokkos) a chance to know something has happened.
4027   */
4028   B->nonzerostate++;
4029   PetscFunctionReturn(PETSC_SUCCESS);
4030 }
4031 
4032 PetscErrorCode MatResetPreallocation_SeqAIJ(Mat A)
4033 {
4034   Mat_SeqAIJ *a;
4035   PetscInt    i;
4036 
4037   PetscFunctionBegin;
4038   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
4039 
4040   /* Check local size. If zero, then return */
4041   if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
4042 
4043   a = (Mat_SeqAIJ *)A->data;
4044   /* if no saved info, we error out */
4045   PetscCheck(a->ipre, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "No saved preallocation info ");
4046 
4047   PetscCheck(a->i && a->j && a->a && a->imax && a->ilen, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Memory info is incomplete, and can not reset preallocation ");
4048 
4049   PetscCall(PetscArraycpy(a->imax, a->ipre, A->rmap->n));
4050   PetscCall(PetscArrayzero(a->ilen, A->rmap->n));
4051   a->i[0] = 0;
4052   for (i = 1; i < A->rmap->n + 1; i++) a->i[i] = a->i[i - 1] + a->imax[i - 1];
4053   A->preallocated     = PETSC_TRUE;
4054   a->nz               = 0;
4055   a->maxnz            = a->i[A->rmap->n];
4056   A->info.nz_unneeded = (double)a->maxnz;
4057   A->was_assembled    = PETSC_FALSE;
4058   A->assembled        = PETSC_FALSE;
4059   PetscFunctionReturn(PETSC_SUCCESS);
4060 }
4061 
4062 /*@
4063    MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in `MATSEQAIJ` format.
4064 
4065    Input Parameters:
4066 +  B - the matrix
4067 .  i - the indices into j for the start of each row (starts with zero)
4068 .  j - the column indices for each row (starts with zero) these must be sorted for each row
4069 -  v - optional values in the matrix
4070 
4071    Level: developer
4072 
4073    Notes:
4074       The `i`,`j`,`v` values are COPIED with this routine; to avoid the copy use `MatCreateSeqAIJWithArrays()`
4075 
4076       This routine may be called multiple times with different nonzero patterns (or the same nonzero pattern). The nonzero
4077       structure will be the union of all the previous nonzero structures.
4078 
4079     Developer Notes:
4080       An optimization could be added to the implementation where it checks if the `i`, and `j` are identical to the current `i` and `j` and
4081       then just copies the `v` values directly with `PetscMemcpy()`.
4082 
4083       This routine could also take a `PetscCopyMode` argument to allow sharing the values instead of always copying them.
4084 
4085 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatSeqAIJSetPreallocation()`, `MatCreateSeqAIJ()`, `MATSEQAIJ`, `MatResetPreallocation()`
4086 @*/
4087 PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B, const PetscInt i[], const PetscInt j[], const PetscScalar v[])
4088 {
4089   PetscFunctionBegin;
4090   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
4091   PetscValidType(B, 1);
4092   PetscTryMethod(B, "MatSeqAIJSetPreallocationCSR_C", (Mat, const PetscInt[], const PetscInt[], const PetscScalar[]), (B, i, j, v));
4093   PetscFunctionReturn(PETSC_SUCCESS);
4094 }
4095 
4096 PetscErrorCode MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B, const PetscInt Ii[], const PetscInt J[], const PetscScalar v[])
4097 {
4098   PetscInt  i;
4099   PetscInt  m, n;
4100   PetscInt  nz;
4101   PetscInt *nnz;
4102 
4103   PetscFunctionBegin;
4104   PetscCheck(Ii[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Ii[0] must be 0 it is %" PetscInt_FMT, Ii[0]);
4105 
4106   PetscCall(PetscLayoutSetUp(B->rmap));
4107   PetscCall(PetscLayoutSetUp(B->cmap));
4108 
4109   PetscCall(MatGetSize(B, &m, &n));
4110   PetscCall(PetscMalloc1(m + 1, &nnz));
4111   for (i = 0; i < m; i++) {
4112     nz = Ii[i + 1] - Ii[i];
4113     PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local row %" PetscInt_FMT " has a negative number of columns %" PetscInt_FMT, i, nz);
4114     nnz[i] = nz;
4115   }
4116   PetscCall(MatSeqAIJSetPreallocation(B, 0, nnz));
4117   PetscCall(PetscFree(nnz));
4118 
4119   for (i = 0; i < m; i++) PetscCall(MatSetValues_SeqAIJ(B, 1, &i, Ii[i + 1] - Ii[i], J + Ii[i], v ? v + Ii[i] : NULL, INSERT_VALUES));
4120 
4121   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
4122   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
4123 
4124   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
4125   PetscFunctionReturn(PETSC_SUCCESS);
4126 }
4127 
4128 /*@
4129    MatSeqAIJKron - Computes `C`, the Kronecker product of `A` and `B`.
4130 
4131    Input Parameters:
4132 +  A - left-hand side matrix
4133 .  B - right-hand side matrix
4134 -  reuse - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
4135 
4136    Output Parameter:
4137 .  C - Kronecker product of `A` and `B`
4138 
4139    Level: intermediate
4140 
4141    Note:
4142       `MAT_REUSE_MATRIX` can only be used when the nonzero structure of the product matrix has not changed from that last call to `MatSeqAIJKron()`.
4143 
4144 .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MATSEQAIJ`, `MATKAIJ`, `MatReuse`
4145 @*/
4146 PetscErrorCode MatSeqAIJKron(Mat A, Mat B, MatReuse reuse, Mat *C)
4147 {
4148   PetscFunctionBegin;
4149   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
4150   PetscValidType(A, 1);
4151   PetscValidHeaderSpecific(B, MAT_CLASSID, 2);
4152   PetscValidType(B, 2);
4153   PetscValidPointer(C, 4);
4154   if (reuse == MAT_REUSE_MATRIX) {
4155     PetscValidHeaderSpecific(*C, MAT_CLASSID, 4);
4156     PetscValidType(*C, 4);
4157   }
4158   PetscTryMethod(A, "MatSeqAIJKron_C", (Mat, Mat, MatReuse, Mat *), (A, B, reuse, C));
4159   PetscFunctionReturn(PETSC_SUCCESS);
4160 }
4161 
4162 PetscErrorCode MatSeqAIJKron_SeqAIJ(Mat A, Mat B, MatReuse reuse, Mat *C)
4163 {
4164   Mat                newmat;
4165   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
4166   Mat_SeqAIJ        *b = (Mat_SeqAIJ *)B->data;
4167   PetscScalar       *v;
4168   const PetscScalar *aa, *ba;
4169   PetscInt          *i, *j, m, n, p, q, nnz = 0, am = A->rmap->n, bm = B->rmap->n, an = A->cmap->n, bn = B->cmap->n;
4170   PetscBool          flg;
4171 
4172   PetscFunctionBegin;
4173   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4174   PetscCheck(A->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
4175   PetscCheck(!B->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4176   PetscCheck(B->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
4177   PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &flg));
4178   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatType %s", ((PetscObject)B)->type_name);
4179   PetscCheck(reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatReuse %d", (int)reuse);
4180   if (reuse == MAT_INITIAL_MATRIX) {
4181     PetscCall(PetscMalloc2(am * bm + 1, &i, a->i[am] * b->i[bm], &j));
4182     PetscCall(MatCreate(PETSC_COMM_SELF, &newmat));
4183     PetscCall(MatSetSizes(newmat, am * bm, an * bn, am * bm, an * bn));
4184     PetscCall(MatSetType(newmat, MATAIJ));
4185     i[0] = 0;
4186     for (m = 0; m < am; ++m) {
4187       for (p = 0; p < bm; ++p) {
4188         i[m * bm + p + 1] = i[m * bm + p] + (a->i[m + 1] - a->i[m]) * (b->i[p + 1] - b->i[p]);
4189         for (n = a->i[m]; n < a->i[m + 1]; ++n) {
4190           for (q = b->i[p]; q < b->i[p + 1]; ++q) j[nnz++] = a->j[n] * bn + b->j[q];
4191         }
4192       }
4193     }
4194     PetscCall(MatSeqAIJSetPreallocationCSR(newmat, i, j, NULL));
4195     *C = newmat;
4196     PetscCall(PetscFree2(i, j));
4197     nnz = 0;
4198   }
4199   PetscCall(MatSeqAIJGetArray(*C, &v));
4200   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
4201   PetscCall(MatSeqAIJGetArrayRead(B, &ba));
4202   for (m = 0; m < am; ++m) {
4203     for (p = 0; p < bm; ++p) {
4204       for (n = a->i[m]; n < a->i[m + 1]; ++n) {
4205         for (q = b->i[p]; q < b->i[p + 1]; ++q) v[nnz++] = aa[n] * ba[q];
4206       }
4207     }
4208   }
4209   PetscCall(MatSeqAIJRestoreArray(*C, &v));
4210   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
4211   PetscCall(MatSeqAIJRestoreArrayRead(B, &ba));
4212   PetscFunctionReturn(PETSC_SUCCESS);
4213 }
4214 
4215 #include <../src/mat/impls/dense/seq/dense.h>
4216 #include <petsc/private/kernels/petscaxpy.h>
4217 
4218 /*
4219     Computes (B'*A')' since computing B*A directly is untenable
4220 
4221                n                       p                          p
4222         [             ]       [             ]         [                 ]
4223       m [      A      ]  *  n [       B     ]   =   m [         C       ]
4224         [             ]       [             ]         [                 ]
4225 
4226 */
4227 PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat A, Mat B, Mat C)
4228 {
4229   Mat_SeqDense      *sub_a = (Mat_SeqDense *)A->data;
4230   Mat_SeqAIJ        *sub_b = (Mat_SeqAIJ *)B->data;
4231   Mat_SeqDense      *sub_c = (Mat_SeqDense *)C->data;
4232   PetscInt           i, j, n, m, q, p;
4233   const PetscInt    *ii, *idx;
4234   const PetscScalar *b, *a, *a_q;
4235   PetscScalar       *c, *c_q;
4236   PetscInt           clda = sub_c->lda;
4237   PetscInt           alda = sub_a->lda;
4238 
4239   PetscFunctionBegin;
4240   m = A->rmap->n;
4241   n = A->cmap->n;
4242   p = B->cmap->n;
4243   a = sub_a->v;
4244   b = sub_b->a;
4245   c = sub_c->v;
4246   if (clda == m) {
4247     PetscCall(PetscArrayzero(c, m * p));
4248   } else {
4249     for (j = 0; j < p; j++)
4250       for (i = 0; i < m; i++) c[j * clda + i] = 0.0;
4251   }
4252   ii  = sub_b->i;
4253   idx = sub_b->j;
4254   for (i = 0; i < n; i++) {
4255     q = ii[i + 1] - ii[i];
4256     while (q-- > 0) {
4257       c_q = c + clda * (*idx);
4258       a_q = a + alda * i;
4259       PetscKernelAXPY(c_q, *b, a_q, m);
4260       idx++;
4261       b++;
4262     }
4263   }
4264   PetscFunctionReturn(PETSC_SUCCESS);
4265 }
4266 
4267 PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat A, Mat B, PetscReal fill, Mat C)
4268 {
4269   PetscInt  m = A->rmap->n, n = B->cmap->n;
4270   PetscBool cisdense;
4271 
4272   PetscFunctionBegin;
4273   PetscCheck(A->cmap->n == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "A->cmap->n %" PetscInt_FMT " != B->rmap->n %" PetscInt_FMT, A->cmap->n, B->rmap->n);
4274   PetscCall(MatSetSizes(C, m, n, m, n));
4275   PetscCall(MatSetBlockSizesFromMats(C, A, B));
4276   PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATSEQDENSECUDA, MATSEQDENSEHIP, ""));
4277   if (!cisdense) PetscCall(MatSetType(C, MATDENSE));
4278   PetscCall(MatSetUp(C));
4279 
4280   C->ops->matmultnumeric = MatMatMultNumeric_SeqDense_SeqAIJ;
4281   PetscFunctionReturn(PETSC_SUCCESS);
4282 }
4283 
4284 /*MC
4285    MATSEQAIJ - MATSEQAIJ = "seqaij" - A matrix type to be used for sequential sparse matrices,
4286    based on compressed sparse row format.
4287 
4288    Options Database Key:
4289 . -mat_type seqaij - sets the matrix type to "seqaij" during a call to MatSetFromOptions()
4290 
4291    Level: beginner
4292 
4293    Notes:
4294     `MatSetValues()` may be called for this matrix type with a `NULL` argument for the numerical values,
4295     in this case the values associated with the rows and columns one passes in are set to zero
4296     in the matrix
4297 
4298     `MatSetOptions`(,`MAT_STRUCTURE_ONLY`,`PETSC_TRUE`) may be called for this matrix type. In this no
4299     space is allocated for the nonzero entries and any entries passed with `MatSetValues()` are ignored
4300 
4301   Developer Note:
4302     It would be nice if all matrix formats supported passing `NULL` in for the numerical values
4303 
4304 .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatSetFromOptions()`, `MatSetType()`, `MatCreate()`, `MatType`, `MATSELL`, `MATSEQSELL`, `MATMPISELL`
4305 M*/
4306 
4307 /*MC
4308    MATAIJ - MATAIJ = "aij" - A matrix type to be used for sparse matrices.
4309 
4310    This matrix type is identical to `MATSEQAIJ` when constructed with a single process communicator,
4311    and `MATMPIAIJ` otherwise.  As a result, for single process communicators,
4312    `MatSeqAIJSetPreallocation()` is supported, and similarly `MatMPIAIJSetPreallocation()` is supported
4313    for communicators controlling multiple processes.  It is recommended that you call both of
4314    the above preallocation routines for simplicity.
4315 
4316    Options Database Key:
4317 . -mat_type aij - sets the matrix type to "aij" during a call to `MatSetFromOptions()`
4318 
4319   Level: beginner
4320 
4321    Note:
4322    Subclasses include `MATAIJCUSPARSE`, `MATAIJPERM`, `MATAIJSELL`, `MATAIJMKL`, `MATAIJCRL`, and also automatically switches over to use inodes when
4323    enough exist.
4324 
4325 .seealso: [](ch_matrices), `Mat`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MATSEQAIJ`, `MATMPIAIJ`, `MATSELL`, `MATSEQSELL`, `MATMPISELL`
4326 M*/
4327 
4328 /*MC
4329    MATAIJCRL - MATAIJCRL = "aijcrl" - A matrix type to be used for sparse matrices.
4330 
4331    Options Database Key:
4332 . -mat_type aijcrl - sets the matrix type to "aijcrl" during a call to `MatSetFromOptions()`
4333 
4334   Level: beginner
4335 
4336    Note:
4337    This matrix type is identical to `MATSEQAIJCRL` when constructed with a single process communicator,
4338    and `MATMPIAIJCRL` otherwise.  As a result, for single process communicators,
4339    `MatSeqAIJSetPreallocation()` is supported, and similarly `MatMPIAIJSetPreallocation()` is supported
4340    for communicators controlling multiple processes.  It is recommended that you call both of
4341    the above preallocation routines for simplicity.
4342 
4343 .seealso: [](ch_matrices), `Mat`, `MatCreateMPIAIJCRL`, `MATSEQAIJCRL`, `MATMPIAIJCRL`, `MATSEQAIJCRL`, `MATMPIAIJCRL`
4344 M*/
4345 
4346 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *);
4347 #if defined(PETSC_HAVE_ELEMENTAL)
4348 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
4349 #endif
4350 #if defined(PETSC_HAVE_SCALAPACK)
4351 PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
4352 #endif
4353 #if defined(PETSC_HAVE_HYPRE)
4354 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType, MatReuse, Mat *);
4355 #endif
4356 
4357 PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat, MatType, MatReuse, Mat *);
4358 PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat, MatType, MatReuse, Mat *);
4359 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_IS_XAIJ(Mat);
4360 
4361 /*@C
4362    MatSeqAIJGetArray - gives read/write access to the array where the data for a `MATSEQAIJ` matrix is stored
4363 
4364    Not Collective
4365 
4366    Input Parameter:
4367 .  mat - a `MATSEQAIJ` matrix
4368 
4369    Output Parameter:
4370 .   array - pointer to the data
4371 
4372    Level: intermediate
4373 
4374    Fortran Note:
4375    `MatSeqAIJGetArray()` Fortran binding is deprecated (since PETSc 3.19), use `MatSeqAIJGetArrayF90()`
4376 
4377 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRestoreArray()`, `MatSeqAIJGetArrayF90()`
4378 @*/
4379 PetscErrorCode MatSeqAIJGetArray(Mat A, PetscScalar **array)
4380 {
4381   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;
4382 
4383   PetscFunctionBegin;
4384   if (aij->ops->getarray) {
4385     PetscCall((*aij->ops->getarray)(A, array));
4386   } else {
4387     *array = aij->a;
4388   }
4389   PetscFunctionReturn(PETSC_SUCCESS);
4390 }
4391 
4392 /*@C
4393    MatSeqAIJRestoreArray - returns access to the array where the data for a `MATSEQAIJ` matrix is stored obtained by `MatSeqAIJGetArray()`
4394 
4395    Not Collective
4396 
4397    Input Parameters:
4398 +  mat - a `MATSEQAIJ` matrix
4399 -  array - pointer to the data
4400 
4401    Level: intermediate
4402 
4403    Fortran Note:
4404    `MatSeqAIJRestoreArray()` Fortran binding is deprecated (since PETSc 3.19), use `MatSeqAIJRestoreArrayF90()`
4405 
4406 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayF90()`
4407 @*/
4408 PetscErrorCode MatSeqAIJRestoreArray(Mat A, PetscScalar **array)
4409 {
4410   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;
4411 
4412   PetscFunctionBegin;
4413   if (aij->ops->restorearray) {
4414     PetscCall((*aij->ops->restorearray)(A, array));
4415   } else {
4416     *array = NULL;
4417   }
4418   PetscCall(MatSeqAIJInvalidateDiagonal(A));
4419   PetscCall(PetscObjectStateIncrease((PetscObject)A));
4420   PetscFunctionReturn(PETSC_SUCCESS);
4421 }
4422 
4423 /*@C
4424    MatSeqAIJGetArrayRead - gives read-only access to the array where the data for a `MATSEQAIJ` matrix is stored
4425 
4426    Not Collective; No Fortran Support
4427 
4428    Input Parameter:
4429 .  mat - a `MATSEQAIJ` matrix
4430 
4431    Output Parameter:
4432 .   array - pointer to the data
4433 
4434    Level: intermediate
4435 
4436 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayRead()`
4437 @*/
4438 PetscErrorCode MatSeqAIJGetArrayRead(Mat A, const PetscScalar **array)
4439 {
4440   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;
4441 
4442   PetscFunctionBegin;
4443   if (aij->ops->getarrayread) {
4444     PetscCall((*aij->ops->getarrayread)(A, array));
4445   } else {
4446     *array = aij->a;
4447   }
4448   PetscFunctionReturn(PETSC_SUCCESS);
4449 }
4450 
4451 /*@C
4452    MatSeqAIJRestoreArrayRead - restore the read-only access array obtained from `MatSeqAIJGetArrayRead()`
4453 
4454    Not Collective; No Fortran Support
4455 
4456    Input Parameter:
4457 .  mat - a `MATSEQAIJ` matrix
4458 
4459    Output Parameter:
4460 .   array - pointer to the data
4461 
4462    Level: intermediate
4463 
4464 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()`
4465 @*/
4466 PetscErrorCode MatSeqAIJRestoreArrayRead(Mat A, const PetscScalar **array)
4467 {
4468   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;
4469 
4470   PetscFunctionBegin;
4471   if (aij->ops->restorearrayread) {
4472     PetscCall((*aij->ops->restorearrayread)(A, array));
4473   } else {
4474     *array = NULL;
4475   }
4476   PetscFunctionReturn(PETSC_SUCCESS);
4477 }
4478 
4479 /*@C
4480    MatSeqAIJGetArrayWrite - gives write-only access to the array where the data for a `MATSEQAIJ` matrix is stored
4481 
4482    Not Collective; No Fortran Support
4483 
4484    Input Parameter:
4485 .  mat - a `MATSEQAIJ` matrix
4486 
4487    Output Parameter:
4488 .   array - pointer to the data
4489 
4490    Level: intermediate
4491 
4492 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayRead()`
4493 @*/
4494 PetscErrorCode MatSeqAIJGetArrayWrite(Mat A, PetscScalar **array)
4495 {
4496   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;
4497 
4498   PetscFunctionBegin;
4499   if (aij->ops->getarraywrite) {
4500     PetscCall((*aij->ops->getarraywrite)(A, array));
4501   } else {
4502     *array = aij->a;
4503   }
4504   PetscCall(MatSeqAIJInvalidateDiagonal(A));
4505   PetscCall(PetscObjectStateIncrease((PetscObject)A));
4506   PetscFunctionReturn(PETSC_SUCCESS);
4507 }
4508 
4509 /*@C
4510    MatSeqAIJRestoreArrayWrite - restore the read-only access array obtained from MatSeqAIJGetArrayRead
4511 
4512    Not Collective; No Fortran Support
4513 
4514    Input Parameter:
4515 .  mat - a MATSEQAIJ matrix
4516 
4517    Output Parameter:
4518 .   array - pointer to the data
4519 
4520    Level: intermediate
4521 
4522 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()`
4523 @*/
4524 PetscErrorCode MatSeqAIJRestoreArrayWrite(Mat A, PetscScalar **array)
4525 {
4526   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;
4527 
4528   PetscFunctionBegin;
4529   if (aij->ops->restorearraywrite) {
4530     PetscCall((*aij->ops->restorearraywrite)(A, array));
4531   } else {
4532     *array = NULL;
4533   }
4534   PetscFunctionReturn(PETSC_SUCCESS);
4535 }
4536 
4537 /*@C
4538    MatSeqAIJGetCSRAndMemType - Get the CSR arrays and the memory type of the `MATSEQAIJ` matrix
4539 
4540    Not Collective; No Fortran Support
4541 
4542    Input Parameter:
4543 .  mat - a matrix of type `MATSEQAIJ` or its subclasses
4544 
4545    Output Parameters:
4546 +  i - row map array of the matrix
4547 .  j - column index array of the matrix
4548 .  a - data array of the matrix
4549 -  memtype - memory type of the arrays
4550 
4551    Level: Developer
4552 
4553   Notes:
4554    Any of the output parameters can be `NULL`, in which case the corresponding value is not returned.
4555    If mat is a device matrix, the arrays are on the device. Otherwise, they are on the host.
4556 
4557    One can call this routine on a preallocated but not assembled matrix to just get the memory of the CSR underneath the matrix.
4558    If the matrix is assembled, the data array `a` is guaranteed to have the latest values of the matrix.
4559 
4560 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()`
4561 @*/
4562 PetscErrorCode MatSeqAIJGetCSRAndMemType(Mat mat, const PetscInt **i, const PetscInt **j, PetscScalar **a, PetscMemType *mtype)
4563 {
4564   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
4565 
4566   PetscFunctionBegin;
4567   PetscCheck(mat->preallocated, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "matrix is not preallocated");
4568   if (aij->ops->getcsrandmemtype) {
4569     PetscCall((*aij->ops->getcsrandmemtype)(mat, i, j, a, mtype));
4570   } else {
4571     if (i) *i = aij->i;
4572     if (j) *j = aij->j;
4573     if (a) *a = aij->a;
4574     if (mtype) *mtype = PETSC_MEMTYPE_HOST;
4575   }
4576   PetscFunctionReturn(PETSC_SUCCESS);
4577 }
4578 
4579 /*@C
4580    MatSeqAIJGetMaxRowNonzeros - returns the maximum number of nonzeros in any row
4581 
4582    Not Collective
4583 
4584    Input Parameter:
4585 .  mat - a `MATSEQAIJ` matrix
4586 
4587    Output Parameter:
4588 .   nz - the maximum number of nonzeros in any row
4589 
4590    Level: intermediate
4591 
4592 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRestoreArray()`, `MatSeqAIJGetArrayF90()`
4593 @*/
4594 PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat A, PetscInt *nz)
4595 {
4596   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;
4597 
4598   PetscFunctionBegin;
4599   *nz = aij->rmax;
4600   PetscFunctionReturn(PETSC_SUCCESS);
4601 }
4602 
4603 static PetscErrorCode MatCOOStructDestroy_SeqAIJ(void *data)
4604 {
4605   MatCOOStruct_SeqAIJ *coo = (MatCOOStruct_SeqAIJ *)data;
4606   PetscFunctionBegin;
4607   PetscCall(PetscFree(coo->perm));
4608   PetscCall(PetscFree(coo->jmap));
4609   PetscCall(PetscFree(coo));
4610   PetscFunctionReturn(PETSC_SUCCESS);
4611 }
4612 
4613 PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat mat, PetscCount coo_n, PetscInt coo_i[], PetscInt coo_j[])
4614 {
4615   MPI_Comm             comm;
4616   PetscInt            *i, *j;
4617   PetscInt             M, N, row;
4618   PetscCount           k, p, q, nneg, nnz, start, end; /* Index the coo array, so use PetscCount as their type */
4619   PetscInt            *Ai;                             /* Change to PetscCount once we use it for row pointers */
4620   PetscInt            *Aj;
4621   PetscScalar         *Aa;
4622   Mat_SeqAIJ          *seqaij = (Mat_SeqAIJ *)(mat->data);
4623   MatType              rtype;
4624   PetscCount          *perm, *jmap;
4625   PetscContainer       container;
4626   MatCOOStruct_SeqAIJ *coo;
4627 
4628   PetscFunctionBegin;
4629   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
4630   PetscCall(MatGetSize(mat, &M, &N));
4631   i = coo_i;
4632   j = coo_j;
4633   PetscCall(PetscMalloc1(coo_n, &perm));
4634   for (k = 0; k < coo_n; k++) { /* Ignore entries with negative row or col indices */
4635     if (j[k] < 0) i[k] = -1;
4636     perm[k] = k;
4637   }
4638 
4639   /* Sort by row */
4640   PetscCall(PetscSortIntWithIntCountArrayPair(coo_n, i, j, perm));
4641   for (k = 0; k < coo_n; k++) {
4642     if (i[k] >= 0) break;
4643   } /* Advance k to the first row with a non-negative index */
4644   nneg = k;
4645   PetscCall(PetscMalloc1(coo_n - nneg + 1, &jmap)); /* +1 to make a CSR-like data structure. jmap[i] originally is the number of repeats for i-th nonzero */
4646   nnz = 0;                                          /* Total number of unique nonzeros to be counted */
4647   jmap++;                                           /* Inc jmap by 1 for convenience */
4648 
4649   PetscCall(PetscCalloc1(M + 1, &Ai));        /* CSR of A */
4650   PetscCall(PetscMalloc1(coo_n - nneg, &Aj)); /* We have at most coo_n-nneg unique nonzeros */
4651 
4652   /* In each row, sort by column, then unique column indices to get row length */
4653   Ai++;  /* Inc by 1 for convenience */
4654   q = 0; /* q-th unique nonzero, with q starting from 0 */
4655   while (k < coo_n) {
4656     row   = i[k];
4657     start = k; /* [start,end) indices for this row */
4658     while (k < coo_n && i[k] == row) k++;
4659     end = k;
4660     PetscCall(PetscSortIntWithCountArray(end - start, j + start, perm + start));
4661     /* Find number of unique col entries in this row */
4662     Aj[q]   = j[start]; /* Log the first nonzero in this row */
4663     jmap[q] = 1;        /* Number of repeats of this nozero entry */
4664     Ai[row] = 1;
4665     nnz++;
4666 
4667     for (p = start + 1; p < end; p++) { /* Scan remaining nonzero in this row */
4668       if (j[p] != j[p - 1]) {           /* Meet a new nonzero */
4669         q++;
4670         jmap[q] = 1;
4671         Aj[q]   = j[p];
4672         Ai[row]++;
4673         nnz++;
4674       } else {
4675         jmap[q]++;
4676       }
4677     }
4678     q++; /* Move to next row and thus next unique nonzero */
4679   }
4680 
4681   Ai--; /* Back to the beginning of Ai[] */
4682   for (k = 0; k < M; k++) Ai[k + 1] += Ai[k];
4683   jmap--; /* Back to the beginning of jmap[] */
4684   jmap[0] = 0;
4685   for (k = 0; k < nnz; k++) jmap[k + 1] += jmap[k];
4686   if (nnz < coo_n - nneg) { /* Realloc with actual number of unique nonzeros */
4687     PetscCount *jmap_new;
4688     PetscInt   *Aj_new;
4689 
4690     PetscCall(PetscMalloc1(nnz + 1, &jmap_new));
4691     PetscCall(PetscArraycpy(jmap_new, jmap, nnz + 1));
4692     PetscCall(PetscFree(jmap));
4693     jmap = jmap_new;
4694 
4695     PetscCall(PetscMalloc1(nnz, &Aj_new));
4696     PetscCall(PetscArraycpy(Aj_new, Aj, nnz));
4697     PetscCall(PetscFree(Aj));
4698     Aj = Aj_new;
4699   }
4700 
4701   if (nneg) { /* Discard heading entries with negative indices in perm[], as we'll access it from index 0 in MatSetValuesCOO */
4702     PetscCount *perm_new;
4703 
4704     PetscCall(PetscMalloc1(coo_n - nneg, &perm_new));
4705     PetscCall(PetscArraycpy(perm_new, perm + nneg, coo_n - nneg));
4706     PetscCall(PetscFree(perm));
4707     perm = perm_new;
4708   }
4709 
4710   PetscCall(MatGetRootType_Private(mat, &rtype));
4711   PetscCall(PetscCalloc1(nnz, &Aa)); /* Zero the matrix */
4712   PetscCall(MatSetSeqAIJWithArrays_private(PETSC_COMM_SELF, M, N, Ai, Aj, Aa, rtype, mat));
4713 
4714   seqaij->singlemalloc = PETSC_FALSE;            /* Ai, Aj and Aa are not allocated in one big malloc */
4715   seqaij->free_a = seqaij->free_ij = PETSC_TRUE; /* Let newmat own Ai, Aj and Aa */
4716 
4717   // Put the COO struct in a container and then attach that to the matrix
4718   PetscCall(PetscMalloc1(1, &coo));
4719   coo->nz   = nnz;
4720   coo->n    = coo_n;
4721   coo->Atot = coo_n - nneg; // Annz is seqaij->nz, so no need to record that again
4722   coo->jmap = jmap;         // of length nnz+1
4723   coo->perm = perm;
4724   PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
4725   PetscCall(PetscContainerSetPointer(container, coo));
4726   PetscCall(PetscContainerSetUserDestroy(container, MatCOOStructDestroy_SeqAIJ));
4727   PetscCall(PetscObjectCompose((PetscObject)mat, "__PETSc_MatCOOStruct_Host", (PetscObject)container));
4728   PetscCall(PetscContainerDestroy(&container));
4729   PetscFunctionReturn(PETSC_SUCCESS);
4730 }
4731 
4732 static PetscErrorCode MatSetValuesCOO_SeqAIJ(Mat A, const PetscScalar v[], InsertMode imode)
4733 {
4734   Mat_SeqAIJ          *aseq = (Mat_SeqAIJ *)A->data;
4735   PetscCount           i, j, Annz = aseq->nz;
4736   PetscCount          *perm, *jmap;
4737   PetscScalar         *Aa;
4738   PetscContainer       container;
4739   MatCOOStruct_SeqAIJ *coo;
4740 
4741   PetscFunctionBegin;
4742   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_MatCOOStruct_Host", (PetscObject *)&container));
4743   PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not found MatCOOStruct on this matrix");
4744   PetscCall(PetscContainerGetPointer(container, (void **)&coo));
4745   perm = coo->perm;
4746   jmap = coo->jmap;
4747   PetscCall(MatSeqAIJGetArray(A, &Aa));
4748   for (i = 0; i < Annz; i++) {
4749     PetscScalar sum = 0.0;
4750     for (j = jmap[i]; j < jmap[i + 1]; j++) sum += v[perm[j]];
4751     Aa[i] = (imode == INSERT_VALUES ? 0.0 : Aa[i]) + sum;
4752   }
4753   PetscCall(MatSeqAIJRestoreArray(A, &Aa));
4754   PetscFunctionReturn(PETSC_SUCCESS);
4755 }
4756 
4757 #if defined(PETSC_HAVE_CUDA)
4758 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCUSPARSE(Mat, MatType, MatReuse, Mat *);
4759 #endif
4760 #if defined(PETSC_HAVE_HIP)
4761 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJHIPSPARSE(Mat, MatType, MatReuse, Mat *);
4762 #endif
4763 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
4764 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJKokkos(Mat, MatType, MatReuse, Mat *);
4765 #endif
4766 
4767 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat B)
4768 {
4769   Mat_SeqAIJ *b;
4770   PetscMPIInt size;
4771 
4772   PetscFunctionBegin;
4773   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
4774   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Comm must be of size 1");
4775 
4776   PetscCall(PetscNew(&b));
4777 
4778   B->data   = (void *)b;
4779   B->ops[0] = MatOps_Values;
4780   if (B->sortedfull) B->ops->setvalues = MatSetValues_SeqAIJ_SortedFull;
4781 
4782   b->row                = NULL;
4783   b->col                = NULL;
4784   b->icol               = NULL;
4785   b->reallocs           = 0;
4786   b->ignorezeroentries  = PETSC_FALSE;
4787   b->roworiented        = PETSC_TRUE;
4788   b->nonew              = 0;
4789   b->diag               = NULL;
4790   b->solve_work         = NULL;
4791   B->spptr              = NULL;
4792   b->saved_values       = NULL;
4793   b->idiag              = NULL;
4794   b->mdiag              = NULL;
4795   b->ssor_work          = NULL;
4796   b->omega              = 1.0;
4797   b->fshift             = 0.0;
4798   b->idiagvalid         = PETSC_FALSE;
4799   b->ibdiagvalid        = PETSC_FALSE;
4800   b->keepnonzeropattern = PETSC_FALSE;
4801 
4802   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ));
4803 #if defined(PETSC_HAVE_MATLAB)
4804   PetscCall(PetscObjectComposeFunction((PetscObject)B, "PetscMatlabEnginePut_C", MatlabEnginePut_SeqAIJ));
4805   PetscCall(PetscObjectComposeFunction((PetscObject)B, "PetscMatlabEngineGet_C", MatlabEngineGet_SeqAIJ));
4806 #endif
4807   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJSetColumnIndices_C", MatSeqAIJSetColumnIndices_SeqAIJ));
4808   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_SeqAIJ));
4809   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_SeqAIJ));
4810   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqsbaij_C", MatConvert_SeqAIJ_SeqSBAIJ));
4811   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqbaij_C", MatConvert_SeqAIJ_SeqBAIJ));
4812   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijperm_C", MatConvert_SeqAIJ_SeqAIJPERM));
4813   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijsell_C", MatConvert_SeqAIJ_SeqAIJSELL));
4814 #if defined(PETSC_HAVE_MKL_SPARSE)
4815   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijmkl_C", MatConvert_SeqAIJ_SeqAIJMKL));
4816 #endif
4817 #if defined(PETSC_HAVE_CUDA)
4818   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijcusparse_C", MatConvert_SeqAIJ_SeqAIJCUSPARSE));
4819   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaijcusparse_seqaij_C", MatProductSetFromOptions_SeqAIJ));
4820   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqaijcusparse_C", MatProductSetFromOptions_SeqAIJ));
4821 #endif
4822 #if defined(PETSC_HAVE_HIP)
4823   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijhipsparse_C", MatConvert_SeqAIJ_SeqAIJHIPSPARSE));
4824   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaijhipsparse_seqaij_C", MatProductSetFromOptions_SeqAIJ));
4825   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqaijhipsparse_C", MatProductSetFromOptions_SeqAIJ));
4826 #endif
4827 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
4828   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijkokkos_C", MatConvert_SeqAIJ_SeqAIJKokkos));
4829 #endif
4830   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijcrl_C", MatConvert_SeqAIJ_SeqAIJCRL));
4831 #if defined(PETSC_HAVE_ELEMENTAL)
4832   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_elemental_C", MatConvert_SeqAIJ_Elemental));
4833 #endif
4834 #if defined(PETSC_HAVE_SCALAPACK)
4835   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_scalapack_C", MatConvert_AIJ_ScaLAPACK));
4836 #endif
4837 #if defined(PETSC_HAVE_HYPRE)
4838   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_hypre_C", MatConvert_AIJ_HYPRE));
4839   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_transpose_seqaij_seqaij_C", MatProductSetFromOptions_Transpose_AIJ_AIJ));
4840 #endif
4841   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqdense_C", MatConvert_SeqAIJ_SeqDense));
4842   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqsell_C", MatConvert_SeqAIJ_SeqSELL));
4843   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_is_C", MatConvert_XAIJ_IS));
4844   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsTranspose_C", MatIsTranspose_SeqAIJ));
4845   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsHermitianTranspose_C", MatIsTranspose_SeqAIJ));
4846   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJSetPreallocation_C", MatSeqAIJSetPreallocation_SeqAIJ));
4847   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatResetPreallocation_C", MatResetPreallocation_SeqAIJ));
4848   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJSetPreallocationCSR_C", MatSeqAIJSetPreallocationCSR_SeqAIJ));
4849   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatReorderForNonzeroDiagonal_C", MatReorderForNonzeroDiagonal_SeqAIJ));
4850   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_is_seqaij_C", MatProductSetFromOptions_IS_XAIJ));
4851   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqdense_seqaij_C", MatProductSetFromOptions_SeqDense_SeqAIJ));
4852   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqaij_C", MatProductSetFromOptions_SeqAIJ));
4853   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJKron_C", MatSeqAIJKron_SeqAIJ));
4854   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSetPreallocationCOO_C", MatSetPreallocationCOO_SeqAIJ));
4855   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSetValuesCOO_C", MatSetValuesCOO_SeqAIJ));
4856   PetscCall(MatCreate_SeqAIJ_Inode(B));
4857   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ));
4858   PetscCall(MatSeqAIJSetTypeFromOptions(B)); /* this allows changing the matrix subtype to say MATSEQAIJPERM */
4859   PetscFunctionReturn(PETSC_SUCCESS);
4860 }
4861 
4862 /*
4863     Given a matrix generated with MatGetFactor() duplicates all the information in A into C
4864 */
4865 PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat C, Mat A, MatDuplicateOption cpvalues, PetscBool mallocmatspace)
4866 {
4867   Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data, *a = (Mat_SeqAIJ *)A->data;
4868   PetscInt    m = A->rmap->n, i;
4869 
4870   PetscFunctionBegin;
4871   PetscCheck(A->assembled || cpvalues == MAT_DO_NOT_COPY_VALUES, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot duplicate unassembled matrix");
4872 
4873   C->factortype = A->factortype;
4874   c->row        = NULL;
4875   c->col        = NULL;
4876   c->icol       = NULL;
4877   c->reallocs   = 0;
4878 
4879   C->assembled = A->assembled;
4880 
4881   if (A->preallocated) {
4882     PetscCall(PetscLayoutReference(A->rmap, &C->rmap));
4883     PetscCall(PetscLayoutReference(A->cmap, &C->cmap));
4884 
4885     if (!A->hash_active) {
4886       PetscCall(PetscMalloc1(m, &c->imax));
4887       PetscCall(PetscMemcpy(c->imax, a->imax, m * sizeof(PetscInt)));
4888       PetscCall(PetscMalloc1(m, &c->ilen));
4889       PetscCall(PetscMemcpy(c->ilen, a->ilen, m * sizeof(PetscInt)));
4890 
4891       /* allocate the matrix space */
4892       if (mallocmatspace) {
4893         PetscCall(PetscMalloc3(a->i[m], &c->a, a->i[m], &c->j, m + 1, &c->i));
4894 
4895         c->singlemalloc = PETSC_TRUE;
4896 
4897         PetscCall(PetscArraycpy(c->i, a->i, m + 1));
4898         if (m > 0) {
4899           PetscCall(PetscArraycpy(c->j, a->j, a->i[m]));
4900           if (cpvalues == MAT_COPY_VALUES) {
4901             const PetscScalar *aa;
4902 
4903             PetscCall(MatSeqAIJGetArrayRead(A, &aa));
4904             PetscCall(PetscArraycpy(c->a, aa, a->i[m]));
4905             PetscCall(MatSeqAIJGetArrayRead(A, &aa));
4906           } else {
4907             PetscCall(PetscArrayzero(c->a, a->i[m]));
4908           }
4909         }
4910       }
4911       C->preallocated = PETSC_TRUE;
4912     } else {
4913       PetscCheck(mallocmatspace, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Cannot malloc matrix memory from a non-preallocated matrix");
4914       PetscCall(MatSetUp(C));
4915     }
4916 
4917     c->ignorezeroentries = a->ignorezeroentries;
4918     c->roworiented       = a->roworiented;
4919     c->nonew             = a->nonew;
4920     if (a->diag) {
4921       PetscCall(PetscMalloc1(m + 1, &c->diag));
4922       PetscCall(PetscMemcpy(c->diag, a->diag, m * sizeof(PetscInt)));
4923     } else c->diag = NULL;
4924 
4925     c->solve_work         = NULL;
4926     c->saved_values       = NULL;
4927     c->idiag              = NULL;
4928     c->ssor_work          = NULL;
4929     c->keepnonzeropattern = a->keepnonzeropattern;
4930     c->free_a             = PETSC_TRUE;
4931     c->free_ij            = PETSC_TRUE;
4932 
4933     c->rmax  = a->rmax;
4934     c->nz    = a->nz;
4935     c->maxnz = a->nz; /* Since we allocate exactly the right amount */
4936 
4937     c->compressedrow.use   = a->compressedrow.use;
4938     c->compressedrow.nrows = a->compressedrow.nrows;
4939     if (a->compressedrow.use) {
4940       i = a->compressedrow.nrows;
4941       PetscCall(PetscMalloc2(i + 1, &c->compressedrow.i, i, &c->compressedrow.rindex));
4942       PetscCall(PetscArraycpy(c->compressedrow.i, a->compressedrow.i, i + 1));
4943       PetscCall(PetscArraycpy(c->compressedrow.rindex, a->compressedrow.rindex, i));
4944     } else {
4945       c->compressedrow.use    = PETSC_FALSE;
4946       c->compressedrow.i      = NULL;
4947       c->compressedrow.rindex = NULL;
4948     }
4949     c->nonzerorowcnt = a->nonzerorowcnt;
4950     C->nonzerostate  = A->nonzerostate;
4951 
4952     PetscCall(MatDuplicate_SeqAIJ_Inode(A, cpvalues, &C));
4953   }
4954   PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist, &((PetscObject)C)->qlist));
4955   PetscFunctionReturn(PETSC_SUCCESS);
4956 }
4957 
4958 PetscErrorCode MatDuplicate_SeqAIJ(Mat A, MatDuplicateOption cpvalues, Mat *B)
4959 {
4960   PetscFunctionBegin;
4961   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
4962   PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n));
4963   if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) PetscCall(MatSetBlockSizesFromMats(*B, A, A));
4964   PetscCall(MatSetType(*B, ((PetscObject)A)->type_name));
4965   PetscCall(MatDuplicateNoCreate_SeqAIJ(*B, A, cpvalues, PETSC_TRUE));
4966   PetscFunctionReturn(PETSC_SUCCESS);
4967 }
4968 
4969 PetscErrorCode MatLoad_SeqAIJ(Mat newMat, PetscViewer viewer)
4970 {
4971   PetscBool isbinary, ishdf5;
4972 
4973   PetscFunctionBegin;
4974   PetscValidHeaderSpecific(newMat, MAT_CLASSID, 1);
4975   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
4976   /* force binary viewer to load .info file if it has not yet done so */
4977   PetscCall(PetscViewerSetUp(viewer));
4978   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
4979   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
4980   if (isbinary) {
4981     PetscCall(MatLoad_SeqAIJ_Binary(newMat, viewer));
4982   } else if (ishdf5) {
4983 #if defined(PETSC_HAVE_HDF5)
4984     PetscCall(MatLoad_AIJ_HDF5(newMat, viewer));
4985 #else
4986     SETERRQ(PetscObjectComm((PetscObject)newMat), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
4987 #endif
4988   } else {
4989     SETERRQ(PetscObjectComm((PetscObject)newMat), PETSC_ERR_SUP, "Viewer type %s not yet supported for reading %s matrices", ((PetscObject)viewer)->type_name, ((PetscObject)newMat)->type_name);
4990   }
4991   PetscFunctionReturn(PETSC_SUCCESS);
4992 }
4993 
4994 PetscErrorCode MatLoad_SeqAIJ_Binary(Mat mat, PetscViewer viewer)
4995 {
4996   Mat_SeqAIJ *a = (Mat_SeqAIJ *)mat->data;
4997   PetscInt    header[4], *rowlens, M, N, nz, sum, rows, cols, i;
4998 
4999   PetscFunctionBegin;
5000   PetscCall(PetscViewerSetUp(viewer));
5001 
5002   /* read in matrix header */
5003   PetscCall(PetscViewerBinaryRead(viewer, header, 4, NULL, PETSC_INT));
5004   PetscCheck(header[0] == MAT_FILE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a matrix object in file");
5005   M  = header[1];
5006   N  = header[2];
5007   nz = header[3];
5008   PetscCheck(M >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Matrix row size (%" PetscInt_FMT ") in file is negative", M);
5009   PetscCheck(N >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Matrix column size (%" PetscInt_FMT ") in file is negative", N);
5010   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Matrix stored in special format on disk, cannot load as SeqAIJ");
5011 
5012   /* set block sizes from the viewer's .info file */
5013   PetscCall(MatLoad_Binary_BlockSizes(mat, viewer));
5014   /* set local and global sizes if not set already */
5015   if (mat->rmap->n < 0) mat->rmap->n = M;
5016   if (mat->cmap->n < 0) mat->cmap->n = N;
5017   if (mat->rmap->N < 0) mat->rmap->N = M;
5018   if (mat->cmap->N < 0) mat->cmap->N = N;
5019   PetscCall(PetscLayoutSetUp(mat->rmap));
5020   PetscCall(PetscLayoutSetUp(mat->cmap));
5021 
5022   /* check if the matrix sizes are correct */
5023   PetscCall(MatGetSize(mat, &rows, &cols));
5024   PetscCheck(M == rows && N == cols, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different sizes (%" PetscInt_FMT ", %" PetscInt_FMT ") than the input matrix (%" PetscInt_FMT ", %" PetscInt_FMT ")", M, N, rows, cols);
5025 
5026   /* read in row lengths */
5027   PetscCall(PetscMalloc1(M, &rowlens));
5028   PetscCall(PetscViewerBinaryRead(viewer, rowlens, M, NULL, PETSC_INT));
5029   /* check if sum(rowlens) is same as nz */
5030   sum = 0;
5031   for (i = 0; i < M; i++) sum += rowlens[i];
5032   PetscCheck(sum == nz, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Inconsistent matrix data in file: nonzeros = %" PetscInt_FMT ", sum-row-lengths = %" PetscInt_FMT, nz, sum);
5033   /* preallocate and check sizes */
5034   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(mat, 0, rowlens));
5035   PetscCall(MatGetSize(mat, &rows, &cols));
5036   PetscCheck(M == rows && N == cols, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%" PetscInt_FMT ", %" PetscInt_FMT ") than the input matrix (%" PetscInt_FMT ", %" PetscInt_FMT ")", M, N, rows, cols);
5037   /* store row lengths */
5038   PetscCall(PetscArraycpy(a->ilen, rowlens, M));
5039   PetscCall(PetscFree(rowlens));
5040 
5041   /* fill in "i" row pointers */
5042   a->i[0] = 0;
5043   for (i = 0; i < M; i++) a->i[i + 1] = a->i[i] + a->ilen[i];
5044   /* read in "j" column indices */
5045   PetscCall(PetscViewerBinaryRead(viewer, a->j, nz, NULL, PETSC_INT));
5046   /* read in "a" nonzero values */
5047   PetscCall(PetscViewerBinaryRead(viewer, a->a, nz, NULL, PETSC_SCALAR));
5048 
5049   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
5050   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
5051   PetscFunctionReturn(PETSC_SUCCESS);
5052 }
5053 
5054 PetscErrorCode MatEqual_SeqAIJ(Mat A, Mat B, PetscBool *flg)
5055 {
5056   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
5057   const PetscScalar *aa, *ba;
5058 #if defined(PETSC_USE_COMPLEX)
5059   PetscInt k;
5060 #endif
5061 
5062   PetscFunctionBegin;
5063   /* If the  matrix dimensions are not equal,or no of nonzeros */
5064   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) || (a->nz != b->nz)) {
5065     *flg = PETSC_FALSE;
5066     PetscFunctionReturn(PETSC_SUCCESS);
5067   }
5068 
5069   /* if the a->i are the same */
5070   PetscCall(PetscArraycmp(a->i, b->i, A->rmap->n + 1, flg));
5071   if (!*flg) PetscFunctionReturn(PETSC_SUCCESS);
5072 
5073   /* if a->j are the same */
5074   PetscCall(PetscArraycmp(a->j, b->j, a->nz, flg));
5075   if (!*flg) PetscFunctionReturn(PETSC_SUCCESS);
5076 
5077   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
5078   PetscCall(MatSeqAIJGetArrayRead(B, &ba));
5079   /* if a->a are the same */
5080 #if defined(PETSC_USE_COMPLEX)
5081   for (k = 0; k < a->nz; k++) {
5082     if (PetscRealPart(aa[k]) != PetscRealPart(ba[k]) || PetscImaginaryPart(aa[k]) != PetscImaginaryPart(ba[k])) {
5083       *flg = PETSC_FALSE;
5084       PetscFunctionReturn(PETSC_SUCCESS);
5085     }
5086   }
5087 #else
5088   PetscCall(PetscArraycmp(aa, ba, a->nz, flg));
5089 #endif
5090   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
5091   PetscCall(MatSeqAIJRestoreArrayRead(B, &ba));
5092   PetscFunctionReturn(PETSC_SUCCESS);
5093 }
5094 
5095 /*@
5096      MatCreateSeqAIJWithArrays - Creates an sequential `MATSEQAIJ` matrix using matrix elements (in CSR format)
5097               provided by the user.
5098 
5099       Collective
5100 
5101    Input Parameters:
5102 +   comm - must be an MPI communicator of size 1
5103 .   m - number of rows
5104 .   n - number of columns
5105 .   i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
5106 .   j - column indices
5107 -   a - matrix values
5108 
5109    Output Parameter:
5110 .   mat - the matrix
5111 
5112    Level: intermediate
5113 
5114    Notes:
5115        The `i`, `j`, and `a` arrays are not copied by this routine, the user must free these arrays
5116     once the matrix is destroyed and not before
5117 
5118        You cannot set new nonzero locations into this matrix, that will generate an error.
5119 
5120        The `i` and `j` indices are 0 based
5121 
5122        The format which is used for the sparse matrix input, is equivalent to a
5123     row-major ordering.. i.e for the following matrix, the input data expected is
5124     as shown
5125 .vb
5126         1 0 0
5127         2 0 3
5128         4 5 6
5129 
5130         i =  {0,1,3,6}  [size = nrow+1  = 3+1]
5131         j =  {0,0,2,0,1,2}  [size = 6]; values must be sorted for each row
5132         v =  {1,2,3,4,5,6}  [size = 6]
5133 .ve
5134 
5135 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MatCreateMPIAIJWithArrays()`, `MatMPIAIJSetPreallocationCSR()`
5136 @*/
5137 PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt i[], PetscInt j[], PetscScalar a[], Mat *mat)
5138 {
5139   PetscInt    ii;
5140   Mat_SeqAIJ *aij;
5141   PetscInt    jj;
5142 
5143   PetscFunctionBegin;
5144   PetscCheck(m <= 0 || i[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "i (row indices) must start with 0");
5145   PetscCall(MatCreate(comm, mat));
5146   PetscCall(MatSetSizes(*mat, m, n, m, n));
5147   /* PetscCall(MatSetBlockSizes(*mat,,)); */
5148   PetscCall(MatSetType(*mat, MATSEQAIJ));
5149   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*mat, MAT_SKIP_ALLOCATION, NULL));
5150   aij = (Mat_SeqAIJ *)(*mat)->data;
5151   PetscCall(PetscMalloc1(m, &aij->imax));
5152   PetscCall(PetscMalloc1(m, &aij->ilen));
5153 
5154   aij->i            = i;
5155   aij->j            = j;
5156   aij->a            = a;
5157   aij->singlemalloc = PETSC_FALSE;
5158   aij->nonew        = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
5159   aij->free_a       = PETSC_FALSE;
5160   aij->free_ij      = PETSC_FALSE;
5161 
5162   for (ii = 0, aij->nonzerorowcnt = 0, aij->rmax = 0; ii < m; ii++) {
5163     aij->ilen[ii] = aij->imax[ii] = i[ii + 1] - i[ii];
5164     if (PetscDefined(USE_DEBUG)) {
5165       PetscCheck(i[ii + 1] - i[ii] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative row length in i (row indices) row = %" PetscInt_FMT " length = %" PetscInt_FMT, ii, i[ii + 1] - i[ii]);
5166       for (jj = i[ii] + 1; jj < i[ii + 1]; jj++) {
5167         PetscCheck(j[jj] >= j[jj - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column entry number %" PetscInt_FMT " (actual column %" PetscInt_FMT ") in row %" PetscInt_FMT " is not sorted", jj - i[ii], j[jj], ii);
5168         PetscCheck(j[jj] != j[jj - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column entry number %" PetscInt_FMT " (actual column %" PetscInt_FMT ") in row %" PetscInt_FMT " is identical to previous entry", jj - i[ii], j[jj], ii);
5169       }
5170     }
5171   }
5172   if (PetscDefined(USE_DEBUG)) {
5173     for (ii = 0; ii < aij->i[m]; ii++) {
5174       PetscCheck(j[ii] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative column index at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]);
5175       PetscCheck(j[ii] <= n - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column index to large at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]);
5176     }
5177   }
5178 
5179   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
5180   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
5181   PetscFunctionReturn(PETSC_SUCCESS);
5182 }
5183 
5184 /*@
5185      MatCreateSeqAIJFromTriple - Creates an sequential `MATSEQAIJ` matrix using matrix elements (in COO format)
5186               provided by the user.
5187 
5188       Collective
5189 
5190    Input Parameters:
5191 +   comm - must be an MPI communicator of size 1
5192 .   m   - number of rows
5193 .   n   - number of columns
5194 .   i   - row indices
5195 .   j   - column indices
5196 .   a   - matrix values
5197 .   nz  - number of nonzeros
5198 -   idx - if the `i` and `j` indices start with 1 use `PETSC_TRUE` otherwise use `PETSC_FALSE`
5199 
5200    Output Parameter:
5201 .   mat - the matrix
5202 
5203    Level: intermediate
5204 
5205    Example:
5206        For the following matrix, the input data expected is as shown (using 0 based indexing)
5207 .vb
5208         1 0 0
5209         2 0 3
5210         4 5 6
5211 
5212         i =  {0,1,1,2,2,2}
5213         j =  {0,0,2,0,1,2}
5214         v =  {1,2,3,4,5,6}
5215 .ve
5216   Note:
5217     Instead of using this function, users should also consider `MatSetPreallocationCOO()` and `MatSetValuesCOO()`, which allow repeated or remote entries,
5218     and are particularly useful in iterative applications.
5219 
5220 .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MatCreateSeqAIJWithArrays()`, `MatMPIAIJSetPreallocationCSR()`, `MatSetValuesCOO()`, `MatSetPreallocationCOO()`
5221 @*/
5222 PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt i[], PetscInt j[], PetscScalar a[], Mat *mat, PetscInt nz, PetscBool idx)
5223 {
5224   PetscInt ii, *nnz, one = 1, row, col;
5225 
5226   PetscFunctionBegin;
5227   PetscCall(PetscCalloc1(m, &nnz));
5228   for (ii = 0; ii < nz; ii++) nnz[i[ii] - !!idx] += 1;
5229   PetscCall(MatCreate(comm, mat));
5230   PetscCall(MatSetSizes(*mat, m, n, m, n));
5231   PetscCall(MatSetType(*mat, MATSEQAIJ));
5232   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*mat, 0, nnz));
5233   for (ii = 0; ii < nz; ii++) {
5234     if (idx) {
5235       row = i[ii] - 1;
5236       col = j[ii] - 1;
5237     } else {
5238       row = i[ii];
5239       col = j[ii];
5240     }
5241     PetscCall(MatSetValues(*mat, one, &row, one, &col, &a[ii], ADD_VALUES));
5242   }
5243   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
5244   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
5245   PetscCall(PetscFree(nnz));
5246   PetscFunctionReturn(PETSC_SUCCESS);
5247 }
5248 
5249 PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat A)
5250 {
5251   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
5252 
5253   PetscFunctionBegin;
5254   a->idiagvalid  = PETSC_FALSE;
5255   a->ibdiagvalid = PETSC_FALSE;
5256 
5257   PetscCall(MatSeqAIJInvalidateDiagonal_Inode(A));
5258   PetscFunctionReturn(PETSC_SUCCESS);
5259 }
5260 
5261 PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm comm, Mat inmat, PetscInt n, MatReuse scall, Mat *outmat)
5262 {
5263   PetscFunctionBegin;
5264   PetscCall(MatCreateMPIMatConcatenateSeqMat_MPIAIJ(comm, inmat, n, scall, outmat));
5265   PetscFunctionReturn(PETSC_SUCCESS);
5266 }
5267 
5268 /*
5269  Permute A into C's *local* index space using rowemb,colemb.
5270  The embedding are supposed to be injections and the above implies that the range of rowemb is a subset
5271  of [0,m), colemb is in [0,n).
5272  If pattern == DIFFERENT_NONZERO_PATTERN, C is preallocated according to A.
5273  */
5274 PetscErrorCode MatSetSeqMat_SeqAIJ(Mat C, IS rowemb, IS colemb, MatStructure pattern, Mat B)
5275 {
5276   /* If making this function public, change the error returned in this function away from _PLIB. */
5277   Mat_SeqAIJ     *Baij;
5278   PetscBool       seqaij;
5279   PetscInt        m, n, *nz, i, j, count;
5280   PetscScalar     v;
5281   const PetscInt *rowindices, *colindices;
5282 
5283   PetscFunctionBegin;
5284   if (!B) PetscFunctionReturn(PETSC_SUCCESS);
5285   /* Check to make sure the target matrix (and embeddings) are compatible with C and each other. */
5286   PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATSEQAIJ, &seqaij));
5287   PetscCheck(seqaij, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is of wrong type");
5288   if (rowemb) {
5289     PetscCall(ISGetLocalSize(rowemb, &m));
5290     PetscCheck(m == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Row IS of size %" PetscInt_FMT " is incompatible with matrix row size %" PetscInt_FMT, m, B->rmap->n);
5291   } else {
5292     PetscCheck(C->rmap->n == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is row-incompatible with the target matrix");
5293   }
5294   if (colemb) {
5295     PetscCall(ISGetLocalSize(colemb, &n));
5296     PetscCheck(n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Diag col IS of size %" PetscInt_FMT " is incompatible with input matrix col size %" PetscInt_FMT, n, B->cmap->n);
5297   } else {
5298     PetscCheck(C->cmap->n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is col-incompatible with the target matrix");
5299   }
5300 
5301   Baij = (Mat_SeqAIJ *)(B->data);
5302   if (pattern == DIFFERENT_NONZERO_PATTERN) {
5303     PetscCall(PetscMalloc1(B->rmap->n, &nz));
5304     for (i = 0; i < B->rmap->n; i++) nz[i] = Baij->i[i + 1] - Baij->i[i];
5305     PetscCall(MatSeqAIJSetPreallocation(C, 0, nz));
5306     PetscCall(PetscFree(nz));
5307   }
5308   if (pattern == SUBSET_NONZERO_PATTERN) PetscCall(MatZeroEntries(C));
5309   count      = 0;
5310   rowindices = NULL;
5311   colindices = NULL;
5312   if (rowemb) PetscCall(ISGetIndices(rowemb, &rowindices));
5313   if (colemb) PetscCall(ISGetIndices(colemb, &colindices));
5314   for (i = 0; i < B->rmap->n; i++) {
5315     PetscInt row;
5316     row = i;
5317     if (rowindices) row = rowindices[i];
5318     for (j = Baij->i[i]; j < Baij->i[i + 1]; j++) {
5319       PetscInt col;
5320       col = Baij->j[count];
5321       if (colindices) col = colindices[col];
5322       v = Baij->a[count];
5323       PetscCall(MatSetValues(C, 1, &row, 1, &col, &v, INSERT_VALUES));
5324       ++count;
5325     }
5326   }
5327   /* FIXME: set C's nonzerostate correctly. */
5328   /* Assembly for C is necessary. */
5329   C->preallocated  = PETSC_TRUE;
5330   C->assembled     = PETSC_TRUE;
5331   C->was_assembled = PETSC_FALSE;
5332   PetscFunctionReturn(PETSC_SUCCESS);
5333 }
5334 
5335 PetscErrorCode MatEliminateZeros_SeqAIJ(Mat A)
5336 {
5337   Mat_SeqAIJ *a  = (Mat_SeqAIJ *)A->data;
5338   MatScalar  *aa = a->a;
5339   PetscInt    m = A->rmap->n, fshift = 0, fshift_prev = 0, i, k;
5340   PetscInt   *ailen = a->ilen, *imax = a->imax, *ai = a->i, *aj = a->j, rmax = 0;
5341 
5342   PetscFunctionBegin;
5343   PetscCheck(A->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot eliminate zeros for unassembled matrix");
5344   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
5345   for (i = 1; i <= m; i++) {
5346     /* move each nonzero entry back by the amount of zero slots (fshift) before it*/
5347     for (k = ai[i - 1]; k < ai[i]; k++) {
5348       if (aa[k] == 0 && aj[k] != i - 1) fshift++;
5349       else {
5350         if (aa[k] == 0 && aj[k] == i - 1) PetscCall(PetscInfo(A, "Keep the diagonal zero at row %" PetscInt_FMT "\n", i - 1));
5351         aa[k - fshift] = aa[k];
5352         aj[k - fshift] = aj[k];
5353       }
5354     }
5355     ai[i - 1] -= fshift_prev; // safe to update ai[i-1] now since it will not be used in the next iteration
5356     fshift_prev = fshift;
5357     /* reset ilen and imax for each row */
5358     ailen[i - 1] = imax[i - 1] = ai[i] - fshift - ai[i - 1];
5359     a->nonzerorowcnt += ((ai[i] - fshift - ai[i - 1]) > 0);
5360     rmax = PetscMax(rmax, ailen[i - 1]);
5361   }
5362   if (m) {
5363     ai[m] -= fshift;
5364     a->nz = ai[m];
5365   }
5366   PetscCall(PetscInfo(A, "Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT "; zeros eliminated: %" PetscInt_FMT "; nonzeros left: %" PetscInt_FMT "\n", m, A->cmap->n, fshift, a->nz));
5367   A->nonzerostate -= fshift;
5368   A->info.nz_unneeded += (PetscReal)fshift;
5369   a->rmax = rmax;
5370   if (a->inode.use && a->inode.checked) PetscCall(MatSeqAIJCheckInode(A));
5371   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5372   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5373   PetscFunctionReturn(PETSC_SUCCESS);
5374 }
5375 
5376 PetscFunctionList MatSeqAIJList = NULL;
5377 
5378 /*@C
5379    MatSeqAIJSetType - Converts a `MATSEQAIJ` matrix to a subtype
5380 
5381    Collective
5382 
5383    Input Parameters:
5384 +  mat      - the matrix object
5385 -  matype   - matrix type
5386 
5387    Options Database Key:
5388 .  -mat_seqaij_type  <method> - for example seqaijcrl
5389 
5390   Level: intermediate
5391 
5392 .seealso: [](ch_matrices), `Mat`, `PCSetType()`, `VecSetType()`, `MatCreate()`, `MatType`, `Mat`
5393 @*/
5394 PetscErrorCode MatSeqAIJSetType(Mat mat, MatType matype)
5395 {
5396   PetscBool sametype;
5397   PetscErrorCode (*r)(Mat, MatType, MatReuse, Mat *);
5398 
5399   PetscFunctionBegin;
5400   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
5401   PetscCall(PetscObjectTypeCompare((PetscObject)mat, matype, &sametype));
5402   if (sametype) PetscFunctionReturn(PETSC_SUCCESS);
5403 
5404   PetscCall(PetscFunctionListFind(MatSeqAIJList, matype, &r));
5405   PetscCheck(r, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Mat type given: %s", matype);
5406   PetscCall((*r)(mat, matype, MAT_INPLACE_MATRIX, &mat));
5407   PetscFunctionReturn(PETSC_SUCCESS);
5408 }
5409 
5410 /*@C
5411   MatSeqAIJRegister -  - Adds a new sub-matrix type for sequential `MATSEQAIJ` matrices
5412 
5413    Not Collective
5414 
5415    Input Parameters:
5416 +  name - name of a new user-defined matrix type, for example `MATSEQAIJCRL`
5417 -  function - routine to convert to subtype
5418 
5419    Level: advanced
5420 
5421    Notes:
5422    `MatSeqAIJRegister()` may be called multiple times to add several user-defined solvers.
5423 
5424    Then, your matrix can be chosen with the procedural interface at runtime via the option
5425 $     -mat_seqaij_type my_mat
5426 
5427 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRegisterAll()`
5428 @*/
5429 PetscErrorCode MatSeqAIJRegister(const char sname[], PetscErrorCode (*function)(Mat, MatType, MatReuse, Mat *))
5430 {
5431   PetscFunctionBegin;
5432   PetscCall(MatInitializePackage());
5433   PetscCall(PetscFunctionListAdd(&MatSeqAIJList, sname, function));
5434   PetscFunctionReturn(PETSC_SUCCESS);
5435 }
5436 
5437 PetscBool MatSeqAIJRegisterAllCalled = PETSC_FALSE;
5438 
5439 /*@C
5440   MatSeqAIJRegisterAll - Registers all of the matrix subtypes of `MATSSEQAIJ`
5441 
5442   Not Collective
5443 
5444   Level: advanced
5445 
5446   Note:
5447   This registers the versions of `MATSEQAIJ` for GPUs
5448 
5449 .seealso: [](ch_matrices), `Mat`, `MatRegisterAll()`, `MatSeqAIJRegister()`
5450 @*/
5451 PetscErrorCode MatSeqAIJRegisterAll(void)
5452 {
5453   PetscFunctionBegin;
5454   if (MatSeqAIJRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
5455   MatSeqAIJRegisterAllCalled = PETSC_TRUE;
5456 
5457   PetscCall(MatSeqAIJRegister(MATSEQAIJCRL, MatConvert_SeqAIJ_SeqAIJCRL));
5458   PetscCall(MatSeqAIJRegister(MATSEQAIJPERM, MatConvert_SeqAIJ_SeqAIJPERM));
5459   PetscCall(MatSeqAIJRegister(MATSEQAIJSELL, MatConvert_SeqAIJ_SeqAIJSELL));
5460 #if defined(PETSC_HAVE_MKL_SPARSE)
5461   PetscCall(MatSeqAIJRegister(MATSEQAIJMKL, MatConvert_SeqAIJ_SeqAIJMKL));
5462 #endif
5463 #if defined(PETSC_HAVE_CUDA)
5464   PetscCall(MatSeqAIJRegister(MATSEQAIJCUSPARSE, MatConvert_SeqAIJ_SeqAIJCUSPARSE));
5465 #endif
5466 #if defined(PETSC_HAVE_HIP)
5467   PetscCall(MatSeqAIJRegister(MATSEQAIJHIPSPARSE, MatConvert_SeqAIJ_SeqAIJHIPSPARSE));
5468 #endif
5469 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
5470   PetscCall(MatSeqAIJRegister(MATSEQAIJKOKKOS, MatConvert_SeqAIJ_SeqAIJKokkos));
5471 #endif
5472 #if defined(PETSC_HAVE_VIENNACL) && defined(PETSC_HAVE_VIENNACL_NO_CUDA)
5473   PetscCall(MatSeqAIJRegister(MATMPIAIJVIENNACL, MatConvert_SeqAIJ_SeqAIJViennaCL));
5474 #endif
5475   PetscFunctionReturn(PETSC_SUCCESS);
5476 }
5477 
5478 /*
5479     Special version for direct calls from Fortran
5480 */
5481 #include <petsc/private/fortranimpl.h>
5482 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5483   #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
5484 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
5485   #define matsetvaluesseqaij_ matsetvaluesseqaij
5486 #endif
5487 
5488 /* Change these macros so can be used in void function */
5489 
5490 /* Change these macros so can be used in void function */
5491 /* Identical to PetscCallVoid, except it assigns to *_ierr */
5492 #undef PetscCall
5493 #define PetscCall(...) \
5494   do { \
5495     PetscErrorCode ierr_msv_mpiaij = __VA_ARGS__; \
5496     if (PetscUnlikely(ierr_msv_mpiaij)) { \
5497       *_ierr = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_msv_mpiaij, PETSC_ERROR_REPEAT, " "); \
5498       return; \
5499     } \
5500   } while (0)
5501 
5502 #undef SETERRQ
5503 #define SETERRQ(comm, ierr, ...) \
5504   do { \
5505     *_ierr = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
5506     return; \
5507   } while (0)
5508 
5509 PETSC_EXTERN void matsetvaluesseqaij_(Mat *AA, PetscInt *mm, const PetscInt im[], PetscInt *nn, const PetscInt in[], const PetscScalar v[], InsertMode *isis, PetscErrorCode *_ierr)
5510 {
5511   Mat         A = *AA;
5512   PetscInt    m = *mm, n = *nn;
5513   InsertMode  is = *isis;
5514   Mat_SeqAIJ *a  = (Mat_SeqAIJ *)A->data;
5515   PetscInt   *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N;
5516   PetscInt   *imax, *ai, *ailen;
5517   PetscInt   *aj, nonew = a->nonew, lastcol = -1;
5518   MatScalar  *ap, value, *aa;
5519   PetscBool   ignorezeroentries = a->ignorezeroentries;
5520   PetscBool   roworiented       = a->roworiented;
5521 
5522   PetscFunctionBegin;
5523   MatCheckPreallocated(A, 1);
5524   imax  = a->imax;
5525   ai    = a->i;
5526   ailen = a->ilen;
5527   aj    = a->j;
5528   aa    = a->a;
5529 
5530   for (k = 0; k < m; k++) { /* loop over added rows */
5531     row = im[k];
5532     if (row < 0) continue;
5533     PetscCheck(row < A->rmap->n, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Row too large");
5534     rp   = aj + ai[row];
5535     ap   = aa + ai[row];
5536     rmax = imax[row];
5537     nrow = ailen[row];
5538     low  = 0;
5539     high = nrow;
5540     for (l = 0; l < n; l++) { /* loop over added columns */
5541       if (in[l] < 0) continue;
5542       PetscCheck(in[l] < A->cmap->n, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Column too large");
5543       col = in[l];
5544       if (roworiented) value = v[l + k * n];
5545       else value = v[k + l * m];
5546 
5547       if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;
5548 
5549       if (col <= lastcol) low = 0;
5550       else high = nrow;
5551       lastcol = col;
5552       while (high - low > 5) {
5553         t = (low + high) / 2;
5554         if (rp[t] > col) high = t;
5555         else low = t;
5556       }
5557       for (i = low; i < high; i++) {
5558         if (rp[i] > col) break;
5559         if (rp[i] == col) {
5560           if (is == ADD_VALUES) ap[i] += value;
5561           else ap[i] = value;
5562           goto noinsert;
5563         }
5564       }
5565       if (value == 0.0 && ignorezeroentries) goto noinsert;
5566       if (nonew == 1) goto noinsert;
5567       PetscCheck(nonew != -1, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero in the matrix");
5568       MatSeqXAIJReallocateAIJ(A, A->rmap->n, 1, nrow, row, col, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar);
5569       N = nrow++ - 1;
5570       a->nz++;
5571       high++;
5572       /* shift up all the later entries in this row */
5573       for (ii = N; ii >= i; ii--) {
5574         rp[ii + 1] = rp[ii];
5575         ap[ii + 1] = ap[ii];
5576       }
5577       rp[i] = col;
5578       ap[i] = value;
5579       A->nonzerostate++;
5580     noinsert:;
5581       low = i + 1;
5582     }
5583     ailen[row] = nrow;
5584   }
5585   PetscFunctionReturnVoid();
5586 }
5587 /* Undefining these here since they were redefined from their original definition above! No
5588  * other PETSc functions should be defined past this point, as it is impossible to recover the
5589  * original definitions */
5590 #undef PetscCall
5591 #undef SETERRQ
5592