Actual source code: aij.c

  1: /*
  2:     Defines the basic matrix operations for the AIJ (compressed row)
  3:   matrix storage format.
  4: */

  6: #include <../src/mat/impls/aij/seq/aij.h>
  7: #include <petscblaslapack.h>
  8: #include <petscbt.h>
  9: #include <petsc/private/kernels/blocktranspose.h>

 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

 19: static PetscErrorCode MatSeqAIJSetTypeFromOptions(Mat A)
 20: {
 21:   PetscBool flg;
 22:   char      type[256];

 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: }

 32: static PetscErrorCode MatGetColumnReductions_SeqAIJ(Mat A, PetscInt type, PetscReal *reductions)
 33: {
 34:   PetscInt    i, m, n;
 35:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

 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");

 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: }

 60: static 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;

 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: }

 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;

 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: }

109: static PetscErrorCode MatFindZeroDiagonals_SeqAIJ(Mat A, IS *zrows)
110: {
111:   PetscInt nrows, *rows;

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: }

120: static 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;

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: }

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;

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: }

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;

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: }

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: }

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;

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++;

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: }

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);

279:   PetscCall(PetscFree(*ia));
280:   PetscCall(PetscFree(*ja));
281:   PetscFunctionReturn(PETSC_SUCCESS);
282: }

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;

297:   PetscFunctionBegin;
298:   *nn = n;
299:   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);

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: }

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: }

335: static 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;

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: }

348: /*
349:     MatSeqAIJSetValuesLocalFast - An optimized version of MatSetValuesLocal() for SeqAIJ matrices with several assumptions

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

356:      This does NOT assume the global column indices are sorted

358: */

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;

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];

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: }

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;

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 = PetscSafePointerPlusOffset(aj, ai[row]);
419:     if (!A->structure_only) ap = PetscSafePointerPlusOffset(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;

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:     noinsert:;
472:     }
473:     ailen[row] = nrow;
474:   }
475:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
476:   PetscFunctionReturn(PETSC_SUCCESS);
477: }

479: static PetscErrorCode MatSetValues_SeqAIJ_SortedFullNoPreallocation(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
480: {
481:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
482:   PetscInt   *rp, k, row;
483:   PetscInt   *ai = a->i;
484:   PetscInt   *aj = a->j;
485:   MatScalar  *aa, *ap;

487:   PetscFunctionBegin;
488:   PetscCheck(!A->was_assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot call on assembled matrix.");
489:   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:   PetscCall(MatSeqAIJGetArray(A, &aa));
492:   for (k = 0; k < m; k++) { /* loop over added rows */
493:     row = im[k];
494:     rp  = aj + ai[row];
495:     ap  = PetscSafePointerPlusOffset(aa, ai[row]);

497:     PetscCall(PetscMemcpy(rp, in, n * sizeof(PetscInt)));
498:     if (!A->structure_only) {
499:       if (v) {
500:         PetscCall(PetscMemcpy(ap, v, n * sizeof(PetscScalar)));
501:         v += n;
502:       } else {
503:         PetscCall(PetscMemzero(ap, n * sizeof(PetscScalar)));
504:       }
505:     }
506:     a->ilen[row]  = n;
507:     a->imax[row]  = n;
508:     a->i[row + 1] = a->i[row] + n;
509:     a->nz += n;
510:   }
511:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
512:   PetscFunctionReturn(PETSC_SUCCESS);
513: }

515: /*@
516:   MatSeqAIJSetTotalPreallocation - Sets an upper bound on the total number of expected nonzeros in the matrix.

518:   Input Parameters:
519: + A       - the `MATSEQAIJ` matrix
520: - nztotal - bound on the number of nonzeros

522:   Level: advanced

524:   Notes:
525:   This can be called if you will be provided the matrix row by row (from row zero) with sorted column indices for each row.
526:   Simply call `MatSetValues()` after this call to provide the matrix entries in the usual manner. This matrix may be used
527:   as always with multiple matrix assemblies.

529: .seealso: [](ch_matrices), `Mat`, `MatSetOption()`, `MAT_SORTED_FULL`, `MatSetValues()`, `MatSeqAIJSetPreallocation()`
530: @*/
531: PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat A, PetscInt nztotal)
532: {
533:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

535:   PetscFunctionBegin;
536:   PetscCall(PetscLayoutSetUp(A->rmap));
537:   PetscCall(PetscLayoutSetUp(A->cmap));
538:   a->maxnz = nztotal;
539:   if (!a->imax) { PetscCall(PetscMalloc1(A->rmap->n, &a->imax)); }
540:   if (!a->ilen) {
541:     PetscCall(PetscMalloc1(A->rmap->n, &a->ilen));
542:   } else {
543:     PetscCall(PetscMemzero(a->ilen, A->rmap->n * sizeof(PetscInt)));
544:   }

546:   /* allocate the matrix space */
547:   PetscCall(PetscShmgetAllocateArray(A->rmap->n + 1, sizeof(PetscInt), (void **)&a->i));
548:   PetscCall(PetscShmgetAllocateArray(nztotal, sizeof(PetscInt), (void **)&a->j));
549:   a->free_ij = PETSC_TRUE;
550:   if (A->structure_only) {
551:     a->free_a = PETSC_FALSE;
552:   } else {
553:     PetscCall(PetscShmgetAllocateArray(nztotal, sizeof(PetscScalar), (void **)&a->a));
554:     a->free_a = PETSC_TRUE;
555:   }
556:   a->i[0]           = 0;
557:   A->ops->setvalues = MatSetValues_SeqAIJ_SortedFullNoPreallocation;
558:   A->preallocated   = PETSC_TRUE;
559:   PetscFunctionReturn(PETSC_SUCCESS);
560: }

562: static PetscErrorCode MatSetValues_SeqAIJ_SortedFull(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
563: {
564:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
565:   PetscInt   *rp, k, row;
566:   PetscInt   *ai = a->i, *ailen = a->ilen;
567:   PetscInt   *aj = a->j;
568:   MatScalar  *aa, *ap;

570:   PetscFunctionBegin;
571:   PetscCall(MatSeqAIJGetArray(A, &aa));
572:   for (k = 0; k < m; k++) { /* loop over added rows */
573:     row = im[k];
574:     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);
575:     rp = aj + ai[row];
576:     ap = aa + ai[row];
577:     if (!A->was_assembled) PetscCall(PetscMemcpy(rp, in, n * sizeof(PetscInt)));
578:     if (!A->structure_only) {
579:       if (v) {
580:         PetscCall(PetscMemcpy(ap, v, n * sizeof(PetscScalar)));
581:         v += n;
582:       } else {
583:         PetscCall(PetscMemzero(ap, n * sizeof(PetscScalar)));
584:       }
585:     }
586:     ailen[row] = n;
587:     a->nz += n;
588:   }
589:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
590:   PetscFunctionReturn(PETSC_SUCCESS);
591: }

593: static PetscErrorCode MatGetValues_SeqAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], PetscScalar v[])
594: {
595:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
596:   PetscInt        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
597:   PetscInt        *ai = a->i, *ailen = a->ilen;
598:   const MatScalar *ap, *aa;

600:   PetscFunctionBegin;
601:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
602:   for (k = 0; k < m; k++) { /* loop over rows */
603:     row = im[k];
604:     if (row < 0) {
605:       v += n;
606:       continue;
607:     } /* negative row */
608:     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);
609:     rp   = PetscSafePointerPlusOffset(aj, ai[row]);
610:     ap   = PetscSafePointerPlusOffset(aa, ai[row]);
611:     nrow = ailen[row];
612:     for (l = 0; l < n; l++) { /* loop over columns */
613:       if (in[l] < 0) {
614:         v++;
615:         continue;
616:       } /* negative column */
617:       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);
618:       col  = in[l];
619:       high = nrow;
620:       low  = 0; /* assume unsorted */
621:       while (high - low > 5) {
622:         t = (low + high) / 2;
623:         if (rp[t] > col) high = t;
624:         else low = t;
625:       }
626:       for (i = low; i < high; i++) {
627:         if (rp[i] > col) break;
628:         if (rp[i] == col) {
629:           *v++ = ap[i];
630:           goto finished;
631:         }
632:       }
633:       *v++ = 0.0;
634:     finished:;
635:     }
636:   }
637:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
638:   PetscFunctionReturn(PETSC_SUCCESS);
639: }

641: static PetscErrorCode MatView_SeqAIJ_Binary(Mat mat, PetscViewer viewer)
642: {
643:   Mat_SeqAIJ        *A = (Mat_SeqAIJ *)mat->data;
644:   const PetscScalar *av;
645:   PetscInt           header[4], M, N, m, nz, i;
646:   PetscInt          *rowlens;

648:   PetscFunctionBegin;
649:   PetscCall(PetscViewerSetUp(viewer));

651:   M  = mat->rmap->N;
652:   N  = mat->cmap->N;
653:   m  = mat->rmap->n;
654:   nz = A->nz;

656:   /* write matrix header */
657:   header[0] = MAT_FILE_CLASSID;
658:   header[1] = M;
659:   header[2] = N;
660:   header[3] = nz;
661:   PetscCall(PetscViewerBinaryWrite(viewer, header, 4, PETSC_INT));

663:   /* fill in and store row lengths */
664:   PetscCall(PetscMalloc1(m, &rowlens));
665:   for (i = 0; i < m; i++) rowlens[i] = A->i[i + 1] - A->i[i];
666:   PetscCall(PetscViewerBinaryWrite(viewer, rowlens, m, PETSC_INT));
667:   PetscCall(PetscFree(rowlens));
668:   /* store column indices */
669:   PetscCall(PetscViewerBinaryWrite(viewer, A->j, nz, PETSC_INT));
670:   /* store nonzero values */
671:   PetscCall(MatSeqAIJGetArrayRead(mat, &av));
672:   PetscCall(PetscViewerBinaryWrite(viewer, av, nz, PETSC_SCALAR));
673:   PetscCall(MatSeqAIJRestoreArrayRead(mat, &av));

675:   /* write block size option to the viewer's .info file */
676:   PetscCall(MatView_Binary_BlockSizes(mat, viewer));
677:   PetscFunctionReturn(PETSC_SUCCESS);
678: }

680: static PetscErrorCode MatView_SeqAIJ_ASCII_structonly(Mat A, PetscViewer viewer)
681: {
682:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
683:   PetscInt    i, k, m = A->rmap->N;

685:   PetscFunctionBegin;
686:   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
687:   for (i = 0; i < m; i++) {
688:     PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
689:     for (k = a->i[i]; k < a->i[i + 1]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ") ", a->j[k]));
690:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
691:   }
692:   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
693:   PetscFunctionReturn(PETSC_SUCCESS);
694: }

696: extern PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat, PetscViewer);

698: static PetscErrorCode MatView_SeqAIJ_ASCII(Mat A, PetscViewer viewer)
699: {
700:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
701:   const PetscScalar *av;
702:   PetscInt           i, j, m = A->rmap->n;
703:   const char        *name;
704:   PetscViewerFormat  format;

706:   PetscFunctionBegin;
707:   if (A->structure_only) {
708:     PetscCall(MatView_SeqAIJ_ASCII_structonly(A, viewer));
709:     PetscFunctionReturn(PETSC_SUCCESS);
710:   }

712:   PetscCall(PetscViewerGetFormat(viewer, &format));
713:   // By petsc's rule, even PETSC_VIEWER_ASCII_INFO_DETAIL doesn't print matrix entries
714:   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS);

716:   /* trigger copy to CPU if needed */
717:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
718:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
719:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
720:     PetscInt nofinalvalue = 0;
721:     if (m && ((a->i[m] == a->i[m - 1]) || (a->j[a->nz - 1] != A->cmap->n - 1))) {
722:       /* Need a dummy value to ensure the dimension of the matrix. */
723:       nofinalvalue = 1;
724:     }
725:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
726:     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Size = %" PetscInt_FMT " %" PetscInt_FMT " \n", m, A->cmap->n));
727:     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Nonzeros = %" PetscInt_FMT " \n", a->nz));
728: #if defined(PETSC_USE_COMPLEX)
729:     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",4);\n", a->nz + nofinalvalue));
730: #else
731:     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",3);\n", a->nz + nofinalvalue));
732: #endif
733:     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = [\n"));

735:     for (i = 0; i < m; i++) {
736:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
737: #if defined(PETSC_USE_COMPLEX)
738:         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])));
739: #else
740:         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n", i + 1, a->j[j] + 1, (double)a->a[j]));
741: #endif
742:       }
743:     }
744:     if (nofinalvalue) {
745: #if defined(PETSC_USE_COMPLEX)
746:       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e %18.16e\n", m, A->cmap->n, 0., 0.));
747: #else
748:       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n", m, A->cmap->n, 0.0));
749: #endif
750:     }
751:     PetscCall(PetscObjectGetName((PetscObject)A, &name));
752:     PetscCall(PetscViewerASCIIPrintf(viewer, "];\n %s = spconvert(zzz);\n", name));
753:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
754:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
755:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
756:     for (i = 0; i < m; i++) {
757:       PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
758:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
759: #if defined(PETSC_USE_COMPLEX)
760:         if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
761:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
762:         } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
763:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)-PetscImaginaryPart(a->a[j])));
764:         } else if (PetscRealPart(a->a[j]) != 0.0) {
765:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j])));
766:         }
767: #else
768:         if (a->a[j] != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j]));
769: #endif
770:       }
771:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
772:     }
773:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
774:   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
775:     PetscInt nzd = 0, fshift = 1, *sptr;
776:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
777:     PetscCall(PetscMalloc1(m + 1, &sptr));
778:     for (i = 0; i < m; i++) {
779:       sptr[i] = nzd + 1;
780:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
781:         if (a->j[j] >= i) {
782: #if defined(PETSC_USE_COMPLEX)
783:           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
784: #else
785:           if (a->a[j] != 0.0) nzd++;
786: #endif
787:         }
788:       }
789:     }
790:     sptr[m] = nzd + 1;
791:     PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT "\n\n", m, nzd));
792:     for (i = 0; i < m + 1; i += 6) {
793:       if (i + 4 < m) {
794:         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]));
795:       } else if (i + 3 < m) {
796:         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]));
797:       } else if (i + 2 < m) {
798:         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2], sptr[i + 3]));
799:       } else if (i + 1 < m) {
800:         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1], sptr[i + 2]));
801:       } else if (i < m) {
802:         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " %" PetscInt_FMT "\n", sptr[i], sptr[i + 1]));
803:       } else {
804:         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "\n", sptr[i]));
805:       }
806:     }
807:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
808:     PetscCall(PetscFree(sptr));
809:     for (i = 0; i < m; i++) {
810:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
811:         if (a->j[j] >= i) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT " ", a->j[j] + fshift));
812:       }
813:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
814:     }
815:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
816:     for (i = 0; i < m; i++) {
817:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
818:         if (a->j[j] >= i) {
819: #if defined(PETSC_USE_COMPLEX)
820:           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])));
821: #else
822:           if (a->a[j] != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " %18.16e ", (double)a->a[j]));
823: #endif
824:         }
825:       }
826:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
827:     }
828:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
829:   } else if (format == PETSC_VIEWER_ASCII_DENSE) {
830:     PetscInt    cnt = 0, jcnt;
831:     PetscScalar value;
832: #if defined(PETSC_USE_COMPLEX)
833:     PetscBool realonly = PETSC_TRUE;

835:     for (i = 0; i < a->i[m]; i++) {
836:       if (PetscImaginaryPart(a->a[i]) != 0.0) {
837:         realonly = PETSC_FALSE;
838:         break;
839:       }
840:     }
841: #endif

843:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
844:     for (i = 0; i < m; i++) {
845:       jcnt = 0;
846:       for (j = 0; j < A->cmap->n; j++) {
847:         if (jcnt < a->i[i + 1] - a->i[i] && j == a->j[cnt]) {
848:           value = a->a[cnt++];
849:           jcnt++;
850:         } else {
851:           value = 0.0;
852:         }
853: #if defined(PETSC_USE_COMPLEX)
854:         if (realonly) {
855:           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)PetscRealPart(value)));
856:         } else {
857:           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e+%7.5e i ", (double)PetscRealPart(value), (double)PetscImaginaryPart(value)));
858:         }
859: #else
860:         PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)value));
861: #endif
862:       }
863:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
864:     }
865:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
866:   } else if (format == PETSC_VIEWER_ASCII_MATRIXMARKET) {
867:     PetscInt fshift = 1;
868:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
869: #if defined(PETSC_USE_COMPLEX)
870:     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate complex general\n"));
871: #else
872:     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate real general\n"));
873: #endif
874:     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", m, A->cmap->n, a->nz));
875:     for (i = 0; i < m; i++) {
876:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
877: #if defined(PETSC_USE_COMPLEX)
878:         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])));
879: #else
880:         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %g\n", i + fshift, a->j[j] + fshift, (double)a->a[j]));
881: #endif
882:       }
883:     }
884:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
885:   } else {
886:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
887:     if (A->factortype) {
888:       for (i = 0; i < m; i++) {
889:         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
890:         /* L part */
891:         for (j = a->i[i]; j < a->i[i + 1]; j++) {
892: #if defined(PETSC_USE_COMPLEX)
893:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
894:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
895:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
896:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)(-PetscImaginaryPart(a->a[j]))));
897:           } else {
898:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j])));
899:           }
900: #else
901:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j]));
902: #endif
903:         }
904:         /* diagonal */
905:         j = a->diag[i];
906: #if defined(PETSC_USE_COMPLEX)
907:         if (PetscImaginaryPart(a->a[j]) > 0.0) {
908:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(1 / a->a[j]), (double)PetscImaginaryPart(1 / a->a[j])));
909:         } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
910:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(1 / a->a[j]), (double)(-PetscImaginaryPart(1 / a->a[j]))));
911:         } else {
912:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(1 / a->a[j])));
913:         }
914: #else
915:         PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)(1 / a->a[j])));
916: #endif

918:         /* U part */
919:         for (j = a->diag[i + 1] + 1; j < a->diag[i]; j++) {
920: #if defined(PETSC_USE_COMPLEX)
921:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
922:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
923:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
924:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)(-PetscImaginaryPart(a->a[j]))));
925:           } else {
926:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j])));
927:           }
928: #else
929:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j]));
930: #endif
931:         }
932:         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
933:       }
934:     } else {
935:       for (i = 0; i < m; i++) {
936:         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
937:         for (j = a->i[i]; j < a->i[i + 1]; j++) {
938: #if defined(PETSC_USE_COMPLEX)
939:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
940:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)PetscImaginaryPart(a->a[j])));
941:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
942:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->j[j], (double)PetscRealPart(a->a[j]), (double)-PetscImaginaryPart(a->a[j])));
943:           } else {
944:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)PetscRealPart(a->a[j])));
945:           }
946: #else
947:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[j], (double)a->a[j]));
948: #endif
949:         }
950:         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
951:       }
952:     }
953:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
954:   }
955:   PetscCall(PetscViewerFlush(viewer));
956:   PetscFunctionReturn(PETSC_SUCCESS);
957: }

959: #include <petscdraw.h>
960: static PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw, void *Aa)
961: {
962:   Mat                A = (Mat)Aa;
963:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
964:   PetscInt           i, j, m = A->rmap->n;
965:   int                color;
966:   PetscReal          xl, yl, xr, yr, x_l, x_r, y_l, y_r;
967:   PetscViewer        viewer;
968:   PetscViewerFormat  format;
969:   const PetscScalar *aa;

971:   PetscFunctionBegin;
972:   PetscCall(PetscObjectQuery((PetscObject)A, "Zoomviewer", (PetscObject *)&viewer));
973:   PetscCall(PetscViewerGetFormat(viewer, &format));
974:   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));

976:   /* loop over matrix elements drawing boxes */
977:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
978:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
979:     PetscDrawCollectiveBegin(draw);
980:     /* Blue for negative, Cyan for zero and  Red for positive */
981:     color = PETSC_DRAW_BLUE;
982:     for (i = 0; i < m; i++) {
983:       y_l = m - i - 1.0;
984:       y_r = y_l + 1.0;
985:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
986:         x_l = a->j[j];
987:         x_r = x_l + 1.0;
988:         if (PetscRealPart(aa[j]) >= 0.) continue;
989:         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
990:       }
991:     }
992:     color = PETSC_DRAW_CYAN;
993:     for (i = 0; i < m; i++) {
994:       y_l = m - i - 1.0;
995:       y_r = y_l + 1.0;
996:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
997:         x_l = a->j[j];
998:         x_r = x_l + 1.0;
999:         if (aa[j] != 0.) continue;
1000:         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1001:       }
1002:     }
1003:     color = PETSC_DRAW_RED;
1004:     for (i = 0; i < m; i++) {
1005:       y_l = m - i - 1.0;
1006:       y_r = y_l + 1.0;
1007:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
1008:         x_l = a->j[j];
1009:         x_r = x_l + 1.0;
1010:         if (PetscRealPart(aa[j]) <= 0.) continue;
1011:         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1012:       }
1013:     }
1014:     PetscDrawCollectiveEnd(draw);
1015:   } else {
1016:     /* use contour shading to indicate magnitude of values */
1017:     /* first determine max of all nonzero values */
1018:     PetscReal minv = 0.0, maxv = 0.0;
1019:     PetscInt  nz = a->nz, count = 0;
1020:     PetscDraw popup;

1022:     for (i = 0; i < nz; i++) {
1023:       if (PetscAbsScalar(aa[i]) > maxv) maxv = PetscAbsScalar(aa[i]);
1024:     }
1025:     if (minv >= maxv) maxv = minv + PETSC_SMALL;
1026:     PetscCall(PetscDrawGetPopup(draw, &popup));
1027:     PetscCall(PetscDrawScalePopup(popup, minv, maxv));

1029:     PetscDrawCollectiveBegin(draw);
1030:     for (i = 0; i < m; i++) {
1031:       y_l = m - i - 1.0;
1032:       y_r = y_l + 1.0;
1033:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
1034:         x_l   = a->j[j];
1035:         x_r   = x_l + 1.0;
1036:         color = PetscDrawRealToColor(PetscAbsScalar(aa[count]), minv, maxv);
1037:         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1038:         count++;
1039:       }
1040:     }
1041:     PetscDrawCollectiveEnd(draw);
1042:   }
1043:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1044:   PetscFunctionReturn(PETSC_SUCCESS);
1045: }

1047: #include <petscdraw.h>
1048: static PetscErrorCode MatView_SeqAIJ_Draw(Mat A, PetscViewer viewer)
1049: {
1050:   PetscDraw draw;
1051:   PetscReal xr, yr, xl, yl, h, w;
1052:   PetscBool isnull;

1054:   PetscFunctionBegin;
1055:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1056:   PetscCall(PetscDrawIsNull(draw, &isnull));
1057:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);

1059:   xr = A->cmap->n;
1060:   yr = A->rmap->n;
1061:   h  = yr / 10.0;
1062:   w  = xr / 10.0;
1063:   xr += w;
1064:   yr += h;
1065:   xl = -w;
1066:   yl = -h;
1067:   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
1068:   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", (PetscObject)viewer));
1069:   PetscCall(PetscDrawZoom(draw, MatView_SeqAIJ_Draw_Zoom, A));
1070:   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", NULL));
1071:   PetscCall(PetscDrawSave(draw));
1072:   PetscFunctionReturn(PETSC_SUCCESS);
1073: }

1075: PetscErrorCode MatView_SeqAIJ(Mat A, PetscViewer viewer)
1076: {
1077:   PetscBool iascii, isbinary, isdraw;

1079:   PetscFunctionBegin;
1080:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1081:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1082:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1083:   if (iascii) PetscCall(MatView_SeqAIJ_ASCII(A, viewer));
1084:   else if (isbinary) PetscCall(MatView_SeqAIJ_Binary(A, viewer));
1085:   else if (isdraw) PetscCall(MatView_SeqAIJ_Draw(A, viewer));
1086:   PetscCall(MatView_SeqAIJ_Inode(A, viewer));
1087:   PetscFunctionReturn(PETSC_SUCCESS);
1088: }

1090: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A, MatAssemblyType mode)
1091: {
1092:   Mat_SeqAIJ *a      = (Mat_SeqAIJ *)A->data;
1093:   PetscInt    fshift = 0, i, *ai = a->i, *aj = a->j, *imax = a->imax;
1094:   PetscInt    m = A->rmap->n, *ip, N, *ailen = a->ilen, rmax = 0, n;
1095:   MatScalar  *aa    = a->a, *ap;
1096:   PetscReal   ratio = 0.6;

1098:   PetscFunctionBegin;
1099:   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
1100:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1101:   if (A->was_assembled && A->ass_nonzerostate == A->nonzerostate) {
1102:     /* we need to respect users asking to use or not the inodes routine in between matrix assemblies */
1103:     PetscCall(MatAssemblyEnd_SeqAIJ_Inode(A, mode));
1104:     PetscFunctionReturn(PETSC_SUCCESS);
1105:   }

1107:   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
1108:   for (i = 1; i < m; i++) {
1109:     /* move each row back by the amount of empty slots (fshift) before it*/
1110:     fshift += imax[i - 1] - ailen[i - 1];
1111:     rmax = PetscMax(rmax, ailen[i]);
1112:     if (fshift) {
1113:       ip = aj + ai[i];
1114:       ap = aa + ai[i];
1115:       N  = ailen[i];
1116:       PetscCall(PetscArraymove(ip - fshift, ip, N));
1117:       if (!A->structure_only) PetscCall(PetscArraymove(ap - fshift, ap, N));
1118:     }
1119:     ai[i] = ai[i - 1] + ailen[i - 1];
1120:   }
1121:   if (m) {
1122:     fshift += imax[m - 1] - ailen[m - 1];
1123:     ai[m] = ai[m - 1] + ailen[m - 1];
1124:   }
1125:   /* reset ilen and imax for each row */
1126:   a->nonzerorowcnt = 0;
1127:   if (A->structure_only) {
1128:     PetscCall(PetscFree(a->imax));
1129:     PetscCall(PetscFree(a->ilen));
1130:   } else { /* !A->structure_only */
1131:     for (i = 0; i < m; i++) {
1132:       ailen[i] = imax[i] = ai[i + 1] - ai[i];
1133:       a->nonzerorowcnt += ((ai[i + 1] - ai[i]) > 0);
1134:     }
1135:   }
1136:   a->nz = ai[m];
1137:   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);
1138:   PetscCall(MatMarkDiagonal_SeqAIJ(A)); // since diagonal info is used a lot, it is helpful to set them up at the end of assembly
1139:   a->diagonaldense = PETSC_TRUE;
1140:   n                = PetscMin(A->rmap->n, A->cmap->n);
1141:   for (i = 0; i < n; i++) {
1142:     if (a->diag[i] >= ai[i + 1]) {
1143:       a->diagonaldense = PETSC_FALSE;
1144:       break;
1145:     }
1146:   }
1147:   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));
1148:   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs));
1149:   PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", rmax));

1151:   A->info.mallocs += a->reallocs;
1152:   a->reallocs         = 0;
1153:   A->info.nz_unneeded = (PetscReal)fshift;
1154:   a->rmax             = rmax;

1156:   if (!A->structure_only) PetscCall(MatCheckCompressedRow(A, a->nonzerorowcnt, &a->compressedrow, a->i, m, ratio));
1157:   PetscCall(MatAssemblyEnd_SeqAIJ_Inode(A, mode));
1158:   PetscFunctionReturn(PETSC_SUCCESS);
1159: }

1161: static PetscErrorCode MatRealPart_SeqAIJ(Mat A)
1162: {
1163:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1164:   PetscInt    i, nz = a->nz;
1165:   MatScalar  *aa;

1167:   PetscFunctionBegin;
1168:   PetscCall(MatSeqAIJGetArray(A, &aa));
1169:   for (i = 0; i < nz; i++) aa[i] = PetscRealPart(aa[i]);
1170:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
1171:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1172:   PetscFunctionReturn(PETSC_SUCCESS);
1173: }

1175: static PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
1176: {
1177:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1178:   PetscInt    i, nz = a->nz;
1179:   MatScalar  *aa;

1181:   PetscFunctionBegin;
1182:   PetscCall(MatSeqAIJGetArray(A, &aa));
1183:   for (i = 0; i < nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1184:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
1185:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1186:   PetscFunctionReturn(PETSC_SUCCESS);
1187: }

1189: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
1190: {
1191:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1192:   MatScalar  *aa;

1194:   PetscFunctionBegin;
1195:   PetscCall(MatSeqAIJGetArrayWrite(A, &aa));
1196:   PetscCall(PetscArrayzero(aa, a->i[A->rmap->n]));
1197:   PetscCall(MatSeqAIJRestoreArrayWrite(A, &aa));
1198:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1199:   PetscFunctionReturn(PETSC_SUCCESS);
1200: }

1202: static PetscErrorCode MatReset_SeqAIJ(Mat A)
1203: {
1204:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

1206:   PetscFunctionBegin;
1207:   if (A->hash_active) {
1208:     A->ops[0] = a->cops;
1209:     PetscCall(PetscHMapIJVDestroy(&a->ht));
1210:     PetscCall(PetscFree(a->dnz));
1211:     A->hash_active = PETSC_FALSE;
1212:   }

1214:   PetscCall(PetscLogObjectState((PetscObject)A, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, A->rmap->n, A->cmap->n, a->nz));
1215:   PetscCall(MatSeqXAIJFreeAIJ(A, &a->a, &a->j, &a->i));
1216:   PetscCall(ISDestroy(&a->row));
1217:   PetscCall(ISDestroy(&a->col));
1218:   PetscCall(PetscFree(a->diag));
1219:   PetscCall(PetscFree(a->ibdiag));
1220:   PetscCall(PetscFree(a->imax));
1221:   PetscCall(PetscFree(a->ilen));
1222:   PetscCall(PetscFree(a->ipre));
1223:   PetscCall(PetscFree3(a->idiag, a->mdiag, a->ssor_work));
1224:   PetscCall(PetscFree(a->solve_work));
1225:   PetscCall(ISDestroy(&a->icol));
1226:   PetscCall(PetscFree(a->saved_values));
1227:   PetscCall(PetscFree2(a->compressedrow.i, a->compressedrow.rindex));
1228:   PetscCall(MatDestroy_SeqAIJ_Inode(A));
1229:   PetscFunctionReturn(PETSC_SUCCESS);
1230: }

1232: static PetscErrorCode MatResetHash_SeqAIJ(Mat A)
1233: {
1234:   PetscFunctionBegin;
1235:   PetscCall(MatReset_SeqAIJ(A));
1236:   PetscCall(MatCreate_SeqAIJ_Inode(A));
1237:   PetscCall(MatSetUp_Seq_Hash(A));
1238:   A->nonzerostate++;
1239:   PetscFunctionReturn(PETSC_SUCCESS);
1240: }

1242: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
1243: {
1244:   PetscFunctionBegin;
1245:   PetscCall(MatReset_SeqAIJ(A));
1246:   PetscCall(PetscFree(A->data));

1248:   /* MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted may allocate this.
1249:      That function is so heavily used (sometimes in an hidden way through multnumeric function pointers)
1250:      that is hard to properly add this data to the MatProduct data. We free it here to avoid
1251:      users reusing the matrix object with different data to incur in obscure segmentation faults
1252:      due to different matrix sizes */
1253:   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc__ab_dense", NULL));

1255:   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
1256:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "PetscMatlabEnginePut_C", NULL));
1257:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "PetscMatlabEngineGet_C", NULL));
1258:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJSetColumnIndices_C", NULL));
1259:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatStoreValues_C", NULL));
1260:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatRetrieveValues_C", NULL));
1261:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqsbaij_C", NULL));
1262:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqbaij_C", NULL));
1263:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijperm_C", NULL));
1264:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijsell_C", NULL));
1265: #if defined(PETSC_HAVE_MKL_SPARSE)
1266:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijmkl_C", NULL));
1267: #endif
1268: #if defined(PETSC_HAVE_CUDA)
1269:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijcusparse_C", NULL));
1270:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijcusparse_seqaij_C", NULL));
1271:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaij_seqaijcusparse_C", NULL));
1272: #endif
1273: #if defined(PETSC_HAVE_HIP)
1274:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijhipsparse_C", NULL));
1275:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijhipsparse_seqaij_C", NULL));
1276:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaij_seqaijhipsparse_C", NULL));
1277: #endif
1278: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1279:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijkokkos_C", NULL));
1280: #endif
1281:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijcrl_C", NULL));
1282: #if defined(PETSC_HAVE_ELEMENTAL)
1283:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_elemental_C", NULL));
1284: #endif
1285: #if defined(PETSC_HAVE_SCALAPACK)
1286:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_scalapack_C", NULL));
1287: #endif
1288: #if defined(PETSC_HAVE_HYPRE)
1289:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_hypre_C", NULL));
1290:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_transpose_seqaij_seqaij_C", NULL));
1291: #endif
1292:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqdense_C", NULL));
1293:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqsell_C", NULL));
1294:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_is_C", NULL));
1295:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatIsTranspose_C", NULL));
1296:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatIsHermitianTranspose_C", NULL));
1297:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", NULL));
1298:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatResetPreallocation_C", NULL));
1299:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatResetHash_C", NULL));
1300:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJSetPreallocationCSR_C", NULL));
1301:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatReorderForNonzeroDiagonal_C", NULL));
1302:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_is_seqaij_C", NULL));
1303:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqdense_seqaij_C", NULL));
1304:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaij_seqaij_C", NULL));
1305:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqAIJKron_C", NULL));
1306:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetPreallocationCOO_C", NULL));
1307:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSetValuesCOO_C", NULL));
1308:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
1309:   /* these calls do not belong here: the subclasses Duplicate/Destroy are wrong */
1310:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaijsell_seqaij_C", NULL));
1311:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaijperm_seqaij_C", NULL));
1312:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaij_seqaijviennacl_C", NULL));
1313:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijviennacl_seqdense_C", NULL));
1314:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_seqaijviennacl_seqaij_C", NULL));
1315:   PetscFunctionReturn(PETSC_SUCCESS);
1316: }

1318: PetscErrorCode MatSetOption_SeqAIJ(Mat A, MatOption op, PetscBool flg)
1319: {
1320:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

1322:   PetscFunctionBegin;
1323:   switch (op) {
1324:   case MAT_ROW_ORIENTED:
1325:     a->roworiented = flg;
1326:     break;
1327:   case MAT_KEEP_NONZERO_PATTERN:
1328:     a->keepnonzeropattern = flg;
1329:     break;
1330:   case MAT_NEW_NONZERO_LOCATIONS:
1331:     a->nonew = (flg ? 0 : 1);
1332:     break;
1333:   case MAT_NEW_NONZERO_LOCATION_ERR:
1334:     a->nonew = (flg ? -1 : 0);
1335:     break;
1336:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1337:     a->nonew = (flg ? -2 : 0);
1338:     break;
1339:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1340:     a->nounused = (flg ? -1 : 0);
1341:     break;
1342:   case MAT_IGNORE_ZERO_ENTRIES:
1343:     a->ignorezeroentries = flg;
1344:     break;
1345:   case MAT_SPD:
1346:   case MAT_SYMMETRIC:
1347:   case MAT_STRUCTURALLY_SYMMETRIC:
1348:   case MAT_HERMITIAN:
1349:   case MAT_SYMMETRY_ETERNAL:
1350:   case MAT_STRUCTURE_ONLY:
1351:   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
1352:   case MAT_SPD_ETERNAL:
1353:     /* if the diagonal matrix is square it inherits some of the properties above */
1354:     break;
1355:   case MAT_FORCE_DIAGONAL_ENTRIES:
1356:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1357:   case MAT_USE_HASH_TABLE:
1358:     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
1359:     break;
1360:   case MAT_USE_INODES:
1361:     PetscCall(MatSetOption_SeqAIJ_Inode(A, MAT_USE_INODES, flg));
1362:     break;
1363:   case MAT_SUBMAT_SINGLEIS:
1364:     A->submat_singleis = flg;
1365:     break;
1366:   case MAT_SORTED_FULL:
1367:     if (flg) A->ops->setvalues = MatSetValues_SeqAIJ_SortedFull;
1368:     else A->ops->setvalues = MatSetValues_SeqAIJ;
1369:     break;
1370:   case MAT_FORM_EXPLICIT_TRANSPOSE:
1371:     A->form_explicit_transpose = flg;
1372:     break;
1373:   default:
1374:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
1375:   }
1376:   PetscFunctionReturn(PETSC_SUCCESS);
1377: }

1379: static PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A, Vec v)
1380: {
1381:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1382:   PetscInt           i, j, n, *ai = a->i, *aj = a->j;
1383:   PetscScalar       *x;
1384:   const PetscScalar *aa;

1386:   PetscFunctionBegin;
1387:   PetscCall(VecGetLocalSize(v, &n));
1388:   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
1389:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
1390:   if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) {
1391:     PetscInt *diag = a->diag;
1392:     PetscCall(VecGetArrayWrite(v, &x));
1393:     for (i = 0; i < n; i++) x[i] = 1.0 / aa[diag[i]];
1394:     PetscCall(VecRestoreArrayWrite(v, &x));
1395:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1396:     PetscFunctionReturn(PETSC_SUCCESS);
1397:   }

1399:   PetscCall(VecGetArrayWrite(v, &x));
1400:   for (i = 0; i < n; i++) {
1401:     x[i] = 0.0;
1402:     for (j = ai[i]; j < ai[i + 1]; j++) {
1403:       if (aj[j] == i) {
1404:         x[i] = aa[j];
1405:         break;
1406:       }
1407:     }
1408:   }
1409:   PetscCall(VecRestoreArrayWrite(v, &x));
1410:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1411:   PetscFunctionReturn(PETSC_SUCCESS);
1412: }

1414: #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>
1415: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A, Vec xx, Vec zz, Vec yy)
1416: {
1417:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1418:   const MatScalar   *aa;
1419:   PetscScalar       *y;
1420:   const PetscScalar *x;
1421:   PetscInt           m = A->rmap->n;
1422: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1423:   const MatScalar  *v;
1424:   PetscScalar       alpha;
1425:   PetscInt          n, i, j;
1426:   const PetscInt   *idx, *ii, *ridx = NULL;
1427:   Mat_CompressedRow cprow    = a->compressedrow;
1428:   PetscBool         usecprow = cprow.use;
1429: #endif

1431:   PetscFunctionBegin;
1432:   if (zz != yy) PetscCall(VecCopy(zz, yy));
1433:   PetscCall(VecGetArrayRead(xx, &x));
1434:   PetscCall(VecGetArray(yy, &y));
1435:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));

1437: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1438:   fortranmulttransposeaddaij_(&m, x, a->i, a->j, aa, y);
1439: #else
1440:   if (usecprow) {
1441:     m    = cprow.nrows;
1442:     ii   = cprow.i;
1443:     ridx = cprow.rindex;
1444:   } else {
1445:     ii = a->i;
1446:   }
1447:   for (i = 0; i < m; i++) {
1448:     idx = a->j + ii[i];
1449:     v   = aa + ii[i];
1450:     n   = ii[i + 1] - ii[i];
1451:     if (usecprow) {
1452:       alpha = x[ridx[i]];
1453:     } else {
1454:       alpha = x[i];
1455:     }
1456:     for (j = 0; j < n; j++) y[idx[j]] += alpha * v[j];
1457:   }
1458: #endif
1459:   PetscCall(PetscLogFlops(2.0 * a->nz));
1460:   PetscCall(VecRestoreArrayRead(xx, &x));
1461:   PetscCall(VecRestoreArray(yy, &y));
1462:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1463:   PetscFunctionReturn(PETSC_SUCCESS);
1464: }

1466: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A, Vec xx, Vec yy)
1467: {
1468:   PetscFunctionBegin;
1469:   PetscCall(VecSet(yy, 0.0));
1470:   PetscCall(MatMultTransposeAdd_SeqAIJ(A, xx, yy, yy));
1471:   PetscFunctionReturn(PETSC_SUCCESS);
1472: }

1474: #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>

1476: PetscErrorCode MatMult_SeqAIJ(Mat A, Vec xx, Vec yy)
1477: {
1478:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1479:   PetscScalar       *y;
1480:   const PetscScalar *x;
1481:   const MatScalar   *a_a;
1482:   PetscInt           m = A->rmap->n;
1483:   const PetscInt    *ii, *ridx = NULL;
1484:   PetscBool          usecprow = a->compressedrow.use;

1486: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1487:   #pragma disjoint(*x, *y, *aa)
1488: #endif

1490:   PetscFunctionBegin;
1491:   if (a->inode.use && a->inode.checked) {
1492:     PetscCall(MatMult_SeqAIJ_Inode(A, xx, yy));
1493:     PetscFunctionReturn(PETSC_SUCCESS);
1494:   }
1495:   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1496:   PetscCall(VecGetArrayRead(xx, &x));
1497:   PetscCall(VecGetArray(yy, &y));
1498:   ii = a->i;
1499:   if (usecprow) { /* use compressed row format */
1500:     PetscCall(PetscArrayzero(y, m));
1501:     m    = a->compressedrow.nrows;
1502:     ii   = a->compressedrow.i;
1503:     ridx = a->compressedrow.rindex;
1504:     PetscPragmaUseOMPKernels(parallel for)
1505:     for (PetscInt i = 0; i < m; i++) {
1506:       PetscInt           n   = ii[i + 1] - ii[i];
1507:       const PetscInt    *aj  = a->j + ii[i];
1508:       const PetscScalar *aa  = a_a + ii[i];
1509:       PetscScalar        sum = 0.0;
1510:       PetscSparseDensePlusDot(sum, x, aa, aj, n);
1511:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1512:       y[*ridx++] = sum;
1513:     }
1514:   } else { /* do not use compressed row format */
1515: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1516:     fortranmultaij_(&m, x, ii, a->j, a_a, y);
1517: #else
1518:     PetscPragmaUseOMPKernels(parallel for)
1519:     for (PetscInt i = 0; i < m; i++) {
1520:       PetscInt           n   = ii[i + 1] - ii[i];
1521:       const PetscInt    *aj  = a->j + ii[i];
1522:       const PetscScalar *aa  = a_a + ii[i];
1523:       PetscScalar        sum = 0.0;
1524:       PetscSparseDensePlusDot(sum, x, aa, aj, n);
1525:       y[i] = sum;
1526:     }
1527: #endif
1528:   }
1529:   PetscCall(PetscLogFlops(2.0 * a->nz - a->nonzerorowcnt));
1530:   PetscCall(VecRestoreArrayRead(xx, &x));
1531:   PetscCall(VecRestoreArray(yy, &y));
1532:   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1533:   PetscFunctionReturn(PETSC_SUCCESS);
1534: }

1536: // HACK!!!!! Used by src/mat/tests/ex170.c
1537: PETSC_EXTERN PetscErrorCode MatMultMax_SeqAIJ(Mat A, Vec xx, Vec yy)
1538: {
1539:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1540:   PetscScalar       *y;
1541:   const PetscScalar *x;
1542:   const MatScalar   *aa, *a_a;
1543:   PetscInt           m = A->rmap->n;
1544:   const PetscInt    *aj, *ii, *ridx   = NULL;
1545:   PetscInt           n, i, nonzerorow = 0;
1546:   PetscScalar        sum;
1547:   PetscBool          usecprow = a->compressedrow.use;

1549: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1550:   #pragma disjoint(*x, *y, *aa)
1551: #endif

1553:   PetscFunctionBegin;
1554:   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1555:   PetscCall(VecGetArrayRead(xx, &x));
1556:   PetscCall(VecGetArray(yy, &y));
1557:   if (usecprow) { /* use compressed row format */
1558:     m    = a->compressedrow.nrows;
1559:     ii   = a->compressedrow.i;
1560:     ridx = a->compressedrow.rindex;
1561:     for (i = 0; i < m; i++) {
1562:       n   = ii[i + 1] - ii[i];
1563:       aj  = a->j + ii[i];
1564:       aa  = a_a + ii[i];
1565:       sum = 0.0;
1566:       nonzerorow += (n > 0);
1567:       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1568:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1569:       y[*ridx++] = sum;
1570:     }
1571:   } else { /* do not use compressed row format */
1572:     ii = a->i;
1573:     for (i = 0; i < m; i++) {
1574:       n   = ii[i + 1] - ii[i];
1575:       aj  = a->j + ii[i];
1576:       aa  = a_a + ii[i];
1577:       sum = 0.0;
1578:       nonzerorow += (n > 0);
1579:       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1580:       y[i] = sum;
1581:     }
1582:   }
1583:   PetscCall(PetscLogFlops(2.0 * a->nz - nonzerorow));
1584:   PetscCall(VecRestoreArrayRead(xx, &x));
1585:   PetscCall(VecRestoreArray(yy, &y));
1586:   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1587:   PetscFunctionReturn(PETSC_SUCCESS);
1588: }

1590: // HACK!!!!! Used by src/mat/tests/ex170.c
1591: PETSC_EXTERN PetscErrorCode MatMultAddMax_SeqAIJ(Mat A, Vec xx, Vec yy, Vec zz)
1592: {
1593:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1594:   PetscScalar       *y, *z;
1595:   const PetscScalar *x;
1596:   const MatScalar   *aa, *a_a;
1597:   PetscInt           m = A->rmap->n, *aj, *ii;
1598:   PetscInt           n, i, *ridx = NULL;
1599:   PetscScalar        sum;
1600:   PetscBool          usecprow = a->compressedrow.use;

1602:   PetscFunctionBegin;
1603:   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1604:   PetscCall(VecGetArrayRead(xx, &x));
1605:   PetscCall(VecGetArrayPair(yy, zz, &y, &z));
1606:   if (usecprow) { /* use compressed row format */
1607:     if (zz != yy) PetscCall(PetscArraycpy(z, y, m));
1608:     m    = a->compressedrow.nrows;
1609:     ii   = a->compressedrow.i;
1610:     ridx = a->compressedrow.rindex;
1611:     for (i = 0; i < m; i++) {
1612:       n   = ii[i + 1] - ii[i];
1613:       aj  = a->j + ii[i];
1614:       aa  = a_a + ii[i];
1615:       sum = y[*ridx];
1616:       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1617:       z[*ridx++] = sum;
1618:     }
1619:   } else { /* do not use compressed row format */
1620:     ii = a->i;
1621:     for (i = 0; i < m; i++) {
1622:       n   = ii[i + 1] - ii[i];
1623:       aj  = a->j + ii[i];
1624:       aa  = a_a + ii[i];
1625:       sum = y[i];
1626:       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1627:       z[i] = sum;
1628:     }
1629:   }
1630:   PetscCall(PetscLogFlops(2.0 * a->nz));
1631:   PetscCall(VecRestoreArrayRead(xx, &x));
1632:   PetscCall(VecRestoreArrayPair(yy, zz, &y, &z));
1633:   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1634:   PetscFunctionReturn(PETSC_SUCCESS);
1635: }

1637: #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h>
1638: PetscErrorCode MatMultAdd_SeqAIJ(Mat A, Vec xx, Vec yy, Vec zz)
1639: {
1640:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1641:   PetscScalar       *y, *z;
1642:   const PetscScalar *x;
1643:   const MatScalar   *a_a;
1644:   const PetscInt    *ii, *ridx = NULL;
1645:   PetscInt           m        = A->rmap->n;
1646:   PetscBool          usecprow = a->compressedrow.use;

1648:   PetscFunctionBegin;
1649:   if (a->inode.use && a->inode.checked) {
1650:     PetscCall(MatMultAdd_SeqAIJ_Inode(A, xx, yy, zz));
1651:     PetscFunctionReturn(PETSC_SUCCESS);
1652:   }
1653:   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1654:   PetscCall(VecGetArrayRead(xx, &x));
1655:   PetscCall(VecGetArrayPair(yy, zz, &y, &z));
1656:   if (usecprow) { /* use compressed row format */
1657:     if (zz != yy) PetscCall(PetscArraycpy(z, y, m));
1658:     m    = a->compressedrow.nrows;
1659:     ii   = a->compressedrow.i;
1660:     ridx = a->compressedrow.rindex;
1661:     for (PetscInt i = 0; i < m; i++) {
1662:       PetscInt           n   = ii[i + 1] - ii[i];
1663:       const PetscInt    *aj  = a->j + ii[i];
1664:       const PetscScalar *aa  = a_a + ii[i];
1665:       PetscScalar        sum = y[*ridx];
1666:       PetscSparseDensePlusDot(sum, x, aa, aj, n);
1667:       z[*ridx++] = sum;
1668:     }
1669:   } else { /* do not use compressed row format */
1670:     ii = a->i;
1671: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1672:     fortranmultaddaij_(&m, x, ii, a->j, a_a, y, z);
1673: #else
1674:     PetscPragmaUseOMPKernels(parallel for)
1675:     for (PetscInt i = 0; i < m; i++) {
1676:       PetscInt           n   = ii[i + 1] - ii[i];
1677:       const PetscInt    *aj  = a->j + ii[i];
1678:       const PetscScalar *aa  = a_a + ii[i];
1679:       PetscScalar        sum = y[i];
1680:       PetscSparseDensePlusDot(sum, x, aa, aj, n);
1681:       z[i] = sum;
1682:     }
1683: #endif
1684:   }
1685:   PetscCall(PetscLogFlops(2.0 * a->nz));
1686:   PetscCall(VecRestoreArrayRead(xx, &x));
1687:   PetscCall(VecRestoreArrayPair(yy, zz, &y, &z));
1688:   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1689:   PetscFunctionReturn(PETSC_SUCCESS);
1690: }

1692: /*
1693:      Adds diagonal pointers to sparse matrix structure.
1694: */
1695: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1696: {
1697:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1698:   PetscInt    i, j, m = A->rmap->n;
1699:   PetscBool   alreadySet = PETSC_TRUE;

1701:   PetscFunctionBegin;
1702:   if (!a->diag) {
1703:     PetscCall(PetscMalloc1(m, &a->diag));
1704:     alreadySet = PETSC_FALSE;
1705:   }
1706:   for (i = 0; i < A->rmap->n; i++) {
1707:     /* If A's diagonal is already correctly set, this fast track enables cheap and repeated MatMarkDiagonal_SeqAIJ() calls */
1708:     if (alreadySet) {
1709:       PetscInt pos = a->diag[i];
1710:       if (pos >= a->i[i] && pos < a->i[i + 1] && a->j[pos] == i) continue;
1711:     }

1713:     a->diag[i] = a->i[i + 1];
1714:     for (j = a->i[i]; j < a->i[i + 1]; j++) {
1715:       if (a->j[j] == i) {
1716:         a->diag[i] = j;
1717:         break;
1718:       }
1719:     }
1720:   }
1721:   PetscFunctionReturn(PETSC_SUCCESS);
1722: }

1724: static PetscErrorCode MatShift_SeqAIJ(Mat A, PetscScalar v)
1725: {
1726:   Mat_SeqAIJ     *a    = (Mat_SeqAIJ *)A->data;
1727:   const PetscInt *diag = (const PetscInt *)a->diag;
1728:   const PetscInt *ii   = (const PetscInt *)a->i;
1729:   PetscInt        i, *mdiag = NULL;
1730:   PetscInt        cnt = 0; /* how many diagonals are missing */

1732:   PetscFunctionBegin;
1733:   if (!A->preallocated || !a->nz) {
1734:     PetscCall(MatSeqAIJSetPreallocation(A, 1, NULL));
1735:     PetscCall(MatShift_Basic(A, v));
1736:     PetscFunctionReturn(PETSC_SUCCESS);
1737:   }

1739:   if (a->diagonaldense) {
1740:     cnt = 0;
1741:   } else {
1742:     PetscCall(PetscCalloc1(A->rmap->n, &mdiag));
1743:     for (i = 0; i < A->rmap->n; i++) {
1744:       if (i < A->cmap->n && diag[i] >= ii[i + 1]) { /* 'out of range' rows never have diagonals */
1745:         cnt++;
1746:         mdiag[i] = 1;
1747:       }
1748:     }
1749:   }
1750:   if (!cnt) {
1751:     PetscCall(MatShift_Basic(A, v));
1752:   } else {
1753:     PetscScalar       *olda = a->a; /* preserve pointers to current matrix nonzeros structure and values */
1754:     PetscInt          *oldj = a->j, *oldi = a->i;
1755:     PetscBool          free_a = a->free_a, free_ij = a->free_ij;
1756:     const PetscScalar *Aa;

1758:     PetscCall(MatSeqAIJGetArrayRead(A, &Aa)); // sync the host
1759:     PetscCall(MatSeqAIJRestoreArrayRead(A, &Aa));

1761:     a->a = NULL;
1762:     a->j = NULL;
1763:     a->i = NULL;
1764:     /* increase the values in imax for each row where a diagonal is being inserted then reallocate the matrix data structures */
1765:     for (i = 0; i < PetscMin(A->rmap->n, A->cmap->n); i++) a->imax[i] += mdiag[i];
1766:     PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(A, 0, a->imax));

1768:     /* copy old values into new matrix data structure */
1769:     for (i = 0; i < A->rmap->n; i++) {
1770:       PetscCall(MatSetValues(A, 1, &i, a->imax[i] - mdiag[i], &oldj[oldi[i]], &olda[oldi[i]], ADD_VALUES));
1771:       if (i < A->cmap->n) PetscCall(MatSetValue(A, i, i, v, ADD_VALUES));
1772:     }
1773:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1774:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1775:     if (free_a) PetscCall(PetscShmgetDeallocateArray((void **)&olda));
1776:     if (free_ij) PetscCall(PetscShmgetDeallocateArray((void **)&oldj));
1777:     if (free_ij) PetscCall(PetscShmgetDeallocateArray((void **)&oldi));
1778:   }
1779:   PetscCall(PetscFree(mdiag));
1780:   a->diagonaldense = PETSC_TRUE;
1781:   PetscFunctionReturn(PETSC_SUCCESS);
1782: }

1784: /*
1785:      Checks for missing diagonals
1786: */
1787: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A, PetscBool *missing, PetscInt *d)
1788: {
1789:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1790:   PetscInt   *diag, *ii = a->i, i;

1792:   PetscFunctionBegin;
1793:   *missing = PETSC_FALSE;
1794:   if (A->rmap->n > 0 && !ii) {
1795:     *missing = PETSC_TRUE;
1796:     if (d) *d = 0;
1797:     PetscCall(PetscInfo(A, "Matrix has no entries therefore is missing diagonal\n"));
1798:   } else {
1799:     PetscInt n;
1800:     n    = PetscMin(A->rmap->n, A->cmap->n);
1801:     diag = a->diag;
1802:     for (i = 0; i < n; i++) {
1803:       if (diag[i] >= ii[i + 1]) {
1804:         *missing = PETSC_TRUE;
1805:         if (d) *d = i;
1806:         PetscCall(PetscInfo(A, "Matrix is missing diagonal number %" PetscInt_FMT "\n", i));
1807:         break;
1808:       }
1809:     }
1810:   }
1811:   PetscFunctionReturn(PETSC_SUCCESS);
1812: }

1814: #include <petscblaslapack.h>
1815: #include <petsc/private/kernels/blockinvert.h>

1817: /*
1818:     Note that values is allocated externally by the PC and then passed into this routine
1819: */
1820: static PetscErrorCode MatInvertVariableBlockDiagonal_SeqAIJ(Mat A, PetscInt nblocks, const PetscInt *bsizes, PetscScalar *diag)
1821: {
1822:   PetscInt        n = A->rmap->n, i, ncnt = 0, *indx, j, bsizemax = 0, *v_pivots;
1823:   PetscBool       allowzeropivot, zeropivotdetected = PETSC_FALSE;
1824:   const PetscReal shift = 0.0;
1825:   PetscInt        ipvt[5];
1826:   PetscCount      flops = 0;
1827:   PetscScalar     work[25], *v_work;

1829:   PetscFunctionBegin;
1830:   allowzeropivot = PetscNot(A->erroriffailure);
1831:   for (i = 0; i < nblocks; i++) ncnt += bsizes[i];
1832:   PetscCheck(ncnt == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Total blocksizes %" PetscInt_FMT " doesn't match number matrix rows %" PetscInt_FMT, ncnt, n);
1833:   for (i = 0; i < nblocks; i++) bsizemax = PetscMax(bsizemax, bsizes[i]);
1834:   PetscCall(PetscMalloc1(bsizemax, &indx));
1835:   if (bsizemax > 7) PetscCall(PetscMalloc2(bsizemax, &v_work, bsizemax, &v_pivots));
1836:   ncnt = 0;
1837:   for (i = 0; i < nblocks; i++) {
1838:     for (j = 0; j < bsizes[i]; j++) indx[j] = ncnt + j;
1839:     PetscCall(MatGetValues(A, bsizes[i], indx, bsizes[i], indx, diag));
1840:     switch (bsizes[i]) {
1841:     case 1:
1842:       *diag = 1.0 / (*diag);
1843:       break;
1844:     case 2:
1845:       PetscCall(PetscKernel_A_gets_inverse_A_2(diag, shift, allowzeropivot, &zeropivotdetected));
1846:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1847:       PetscCall(PetscKernel_A_gets_transpose_A_2(diag));
1848:       break;
1849:     case 3:
1850:       PetscCall(PetscKernel_A_gets_inverse_A_3(diag, shift, allowzeropivot, &zeropivotdetected));
1851:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1852:       PetscCall(PetscKernel_A_gets_transpose_A_3(diag));
1853:       break;
1854:     case 4:
1855:       PetscCall(PetscKernel_A_gets_inverse_A_4(diag, shift, allowzeropivot, &zeropivotdetected));
1856:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1857:       PetscCall(PetscKernel_A_gets_transpose_A_4(diag));
1858:       break;
1859:     case 5:
1860:       PetscCall(PetscKernel_A_gets_inverse_A_5(diag, ipvt, work, shift, allowzeropivot, &zeropivotdetected));
1861:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1862:       PetscCall(PetscKernel_A_gets_transpose_A_5(diag));
1863:       break;
1864:     case 6:
1865:       PetscCall(PetscKernel_A_gets_inverse_A_6(diag, shift, allowzeropivot, &zeropivotdetected));
1866:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1867:       PetscCall(PetscKernel_A_gets_transpose_A_6(diag));
1868:       break;
1869:     case 7:
1870:       PetscCall(PetscKernel_A_gets_inverse_A_7(diag, shift, allowzeropivot, &zeropivotdetected));
1871:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1872:       PetscCall(PetscKernel_A_gets_transpose_A_7(diag));
1873:       break;
1874:     default:
1875:       PetscCall(PetscKernel_A_gets_inverse_A(bsizes[i], diag, v_pivots, v_work, allowzeropivot, &zeropivotdetected));
1876:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1877:       PetscCall(PetscKernel_A_gets_transpose_A_N(diag, bsizes[i]));
1878:     }
1879:     ncnt += bsizes[i];
1880:     diag += bsizes[i] * bsizes[i];
1881:     flops += 2 * PetscPowInt64(bsizes[i], 3) / 3;
1882:   }
1883:   PetscCall(PetscLogFlops(flops));
1884:   if (bsizemax > 7) PetscCall(PetscFree2(v_work, v_pivots));
1885:   PetscCall(PetscFree(indx));
1886:   PetscFunctionReturn(PETSC_SUCCESS);
1887: }

1889: /*
1890:    Negative shift indicates do not generate an error if there is a zero diagonal, just invert it anyways
1891: */
1892: static PetscErrorCode MatInvertDiagonal_SeqAIJ(Mat A, PetscScalar omega, PetscScalar fshift)
1893: {
1894:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
1895:   PetscInt         i, *diag, m = A->rmap->n;
1896:   const MatScalar *v;
1897:   PetscScalar     *idiag, *mdiag;

1899:   PetscFunctionBegin;
1900:   if (a->idiagvalid) PetscFunctionReturn(PETSC_SUCCESS);
1901:   PetscCall(MatMarkDiagonal_SeqAIJ(A));
1902:   diag = a->diag;
1903:   if (!a->idiag) { PetscCall(PetscMalloc3(m, &a->idiag, m, &a->mdiag, m, &a->ssor_work)); }

1905:   mdiag = a->mdiag;
1906:   idiag = a->idiag;
1907:   PetscCall(MatSeqAIJGetArrayRead(A, &v));
1908:   if (omega == 1.0 && PetscRealPart(fshift) <= 0.0) {
1909:     for (i = 0; i < m; i++) {
1910:       mdiag[i] = v[diag[i]];
1911:       if (!PetscAbsScalar(mdiag[i])) { /* zero diagonal */
1912:         if (PetscRealPart(fshift)) {
1913:           PetscCall(PetscInfo(A, "Zero diagonal on row %" PetscInt_FMT "\n", i));
1914:           A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1915:           A->factorerror_zeropivot_value = 0.0;
1916:           A->factorerror_zeropivot_row   = i;
1917:         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Zero diagonal on row %" PetscInt_FMT, i);
1918:       }
1919:       idiag[i] = 1.0 / v[diag[i]];
1920:     }
1921:     PetscCall(PetscLogFlops(m));
1922:   } else {
1923:     for (i = 0; i < m; i++) {
1924:       mdiag[i] = v[diag[i]];
1925:       idiag[i] = omega / (fshift + v[diag[i]]);
1926:     }
1927:     PetscCall(PetscLogFlops(2.0 * m));
1928:   }
1929:   a->idiagvalid = PETSC_TRUE;
1930:   PetscCall(MatSeqAIJRestoreArrayRead(A, &v));
1931:   PetscFunctionReturn(PETSC_SUCCESS);
1932: }

1934: PetscErrorCode MatSOR_SeqAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
1935: {
1936:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1937:   PetscScalar       *x, d, sum, *t, scale;
1938:   const MatScalar   *v, *idiag = NULL, *mdiag, *aa;
1939:   const PetscScalar *b, *bs, *xb, *ts;
1940:   PetscInt           n, m = A->rmap->n, i;
1941:   const PetscInt    *idx, *diag;

1943:   PetscFunctionBegin;
1944:   if (a->inode.use && a->inode.checked && omega == 1.0 && fshift == 0.0) {
1945:     PetscCall(MatSOR_SeqAIJ_Inode(A, bb, omega, flag, fshift, its, lits, xx));
1946:     PetscFunctionReturn(PETSC_SUCCESS);
1947:   }
1948:   its = its * lits;

1950:   if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
1951:   if (!a->idiagvalid) PetscCall(MatInvertDiagonal_SeqAIJ(A, omega, fshift));
1952:   a->fshift = fshift;
1953:   a->omega  = omega;

1955:   diag  = a->diag;
1956:   t     = a->ssor_work;
1957:   idiag = a->idiag;
1958:   mdiag = a->mdiag;

1960:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
1961:   PetscCall(VecGetArray(xx, &x));
1962:   PetscCall(VecGetArrayRead(bb, &b));
1963:   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1964:   if (flag == SOR_APPLY_UPPER) {
1965:     /* apply (U + D/omega) to the vector */
1966:     bs = b;
1967:     for (i = 0; i < m; i++) {
1968:       d   = fshift + mdiag[i];
1969:       n   = a->i[i + 1] - diag[i] - 1;
1970:       idx = a->j + diag[i] + 1;
1971:       v   = aa + diag[i] + 1;
1972:       sum = b[i] * d / omega;
1973:       PetscSparseDensePlusDot(sum, bs, v, idx, n);
1974:       x[i] = sum;
1975:     }
1976:     PetscCall(VecRestoreArray(xx, &x));
1977:     PetscCall(VecRestoreArrayRead(bb, &b));
1978:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1979:     PetscCall(PetscLogFlops(a->nz));
1980:     PetscFunctionReturn(PETSC_SUCCESS);
1981:   }

1983:   PetscCheck(flag != SOR_APPLY_LOWER, PETSC_COMM_SELF, PETSC_ERR_SUP, "SOR_APPLY_LOWER is not implemented");
1984:   if (flag & SOR_EISENSTAT) {
1985:     /* Let  A = L + U + D; where L is lower triangular,
1986:     U is upper triangular, E = D/omega; This routine applies

1988:             (L + E)^{-1} A (U + E)^{-1}

1990:     to a vector efficiently using Eisenstat's trick.
1991:     */
1992:     scale = (2.0 / omega) - 1.0;

1994:     /*  x = (E + U)^{-1} b */
1995:     for (i = m - 1; i >= 0; i--) {
1996:       n   = a->i[i + 1] - diag[i] - 1;
1997:       idx = a->j + diag[i] + 1;
1998:       v   = aa + diag[i] + 1;
1999:       sum = b[i];
2000:       PetscSparseDenseMinusDot(sum, x, v, idx, n);
2001:       x[i] = sum * idiag[i];
2002:     }

2004:     /*  t = b - (2*E - D)x */
2005:     v = aa;
2006:     for (i = 0; i < m; i++) t[i] = b[i] - scale * (v[*diag++]) * x[i];

2008:     /*  t = (E + L)^{-1}t */
2009:     ts   = t;
2010:     diag = a->diag;
2011:     for (i = 0; i < m; i++) {
2012:       n   = diag[i] - a->i[i];
2013:       idx = a->j + a->i[i];
2014:       v   = aa + a->i[i];
2015:       sum = t[i];
2016:       PetscSparseDenseMinusDot(sum, ts, v, idx, n);
2017:       t[i] = sum * idiag[i];
2018:       /*  x = x + t */
2019:       x[i] += t[i];
2020:     }

2022:     PetscCall(PetscLogFlops(6.0 * m - 1 + 2.0 * a->nz));
2023:     PetscCall(VecRestoreArray(xx, &x));
2024:     PetscCall(VecRestoreArrayRead(bb, &b));
2025:     PetscFunctionReturn(PETSC_SUCCESS);
2026:   }
2027:   if (flag & SOR_ZERO_INITIAL_GUESS) {
2028:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
2029:       for (i = 0; i < m; i++) {
2030:         n   = diag[i] - a->i[i];
2031:         idx = a->j + a->i[i];
2032:         v   = aa + a->i[i];
2033:         sum = b[i];
2034:         PetscSparseDenseMinusDot(sum, x, v, idx, n);
2035:         t[i] = sum;
2036:         x[i] = sum * idiag[i];
2037:       }
2038:       xb = t;
2039:       PetscCall(PetscLogFlops(a->nz));
2040:     } else xb = b;
2041:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
2042:       for (i = m - 1; i >= 0; i--) {
2043:         n   = a->i[i + 1] - diag[i] - 1;
2044:         idx = a->j + diag[i] + 1;
2045:         v   = aa + diag[i] + 1;
2046:         sum = xb[i];
2047:         PetscSparseDenseMinusDot(sum, x, v, idx, n);
2048:         if (xb == b) {
2049:           x[i] = sum * idiag[i];
2050:         } else {
2051:           x[i] = (1 - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
2052:         }
2053:       }
2054:       PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
2055:     }
2056:     its--;
2057:   }
2058:   while (its--) {
2059:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
2060:       for (i = 0; i < m; i++) {
2061:         /* lower */
2062:         n   = diag[i] - a->i[i];
2063:         idx = a->j + a->i[i];
2064:         v   = aa + a->i[i];
2065:         sum = b[i];
2066:         PetscSparseDenseMinusDot(sum, x, v, idx, n);
2067:         t[i] = sum; /* save application of the lower-triangular part */
2068:         /* upper */
2069:         n   = a->i[i + 1] - diag[i] - 1;
2070:         idx = a->j + diag[i] + 1;
2071:         v   = aa + diag[i] + 1;
2072:         PetscSparseDenseMinusDot(sum, x, v, idx, n);
2073:         x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
2074:       }
2075:       xb = t;
2076:       PetscCall(PetscLogFlops(2.0 * a->nz));
2077:     } else xb = b;
2078:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
2079:       for (i = m - 1; i >= 0; i--) {
2080:         sum = xb[i];
2081:         if (xb == b) {
2082:           /* whole matrix (no checkpointing available) */
2083:           n   = a->i[i + 1] - a->i[i];
2084:           idx = a->j + a->i[i];
2085:           v   = aa + a->i[i];
2086:           PetscSparseDenseMinusDot(sum, x, v, idx, n);
2087:           x[i] = (1. - omega) * x[i] + (sum + mdiag[i] * x[i]) * idiag[i];
2088:         } else { /* lower-triangular part has been saved, so only apply upper-triangular */
2089:           n   = a->i[i + 1] - diag[i] - 1;
2090:           idx = a->j + diag[i] + 1;
2091:           v   = aa + diag[i] + 1;
2092:           PetscSparseDenseMinusDot(sum, x, v, idx, n);
2093:           x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
2094:         }
2095:       }
2096:       if (xb == b) {
2097:         PetscCall(PetscLogFlops(2.0 * a->nz));
2098:       } else {
2099:         PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
2100:       }
2101:     }
2102:   }
2103:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2104:   PetscCall(VecRestoreArray(xx, &x));
2105:   PetscCall(VecRestoreArrayRead(bb, &b));
2106:   PetscFunctionReturn(PETSC_SUCCESS);
2107: }

2109: static PetscErrorCode MatGetInfo_SeqAIJ(Mat A, MatInfoType flag, MatInfo *info)
2110: {
2111:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

2113:   PetscFunctionBegin;
2114:   info->block_size   = 1.0;
2115:   info->nz_allocated = a->maxnz;
2116:   info->nz_used      = a->nz;
2117:   info->nz_unneeded  = (a->maxnz - a->nz);
2118:   info->assemblies   = A->num_ass;
2119:   info->mallocs      = A->info.mallocs;
2120:   info->memory       = 0; /* REVIEW ME */
2121:   if (A->factortype) {
2122:     info->fill_ratio_given  = A->info.fill_ratio_given;
2123:     info->fill_ratio_needed = A->info.fill_ratio_needed;
2124:     info->factor_mallocs    = A->info.factor_mallocs;
2125:   } else {
2126:     info->fill_ratio_given  = 0;
2127:     info->fill_ratio_needed = 0;
2128:     info->factor_mallocs    = 0;
2129:   }
2130:   PetscFunctionReturn(PETSC_SUCCESS);
2131: }

2133: static PetscErrorCode MatZeroRows_SeqAIJ(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
2134: {
2135:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2136:   PetscInt           i, m = A->rmap->n - 1;
2137:   const PetscScalar *xx;
2138:   PetscScalar       *bb, *aa;
2139:   PetscInt           d = 0;

2141:   PetscFunctionBegin;
2142:   if (x && b) {
2143:     PetscCall(VecGetArrayRead(x, &xx));
2144:     PetscCall(VecGetArray(b, &bb));
2145:     for (i = 0; i < N; i++) {
2146:       PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2147:       if (rows[i] >= A->cmap->n) continue;
2148:       bb[rows[i]] = diag * xx[rows[i]];
2149:     }
2150:     PetscCall(VecRestoreArrayRead(x, &xx));
2151:     PetscCall(VecRestoreArray(b, &bb));
2152:   }

2154:   PetscCall(MatSeqAIJGetArray(A, &aa));
2155:   if (a->keepnonzeropattern) {
2156:     for (i = 0; i < N; i++) {
2157:       PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2158:       PetscCall(PetscArrayzero(&aa[a->i[rows[i]]], a->ilen[rows[i]]));
2159:     }
2160:     if (diag != 0.0) {
2161:       for (i = 0; i < N; i++) {
2162:         d = rows[i];
2163:         if (rows[i] >= A->cmap->n) continue;
2164:         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);
2165:       }
2166:       for (i = 0; i < N; i++) {
2167:         if (rows[i] >= A->cmap->n) continue;
2168:         aa[a->diag[rows[i]]] = diag;
2169:       }
2170:     }
2171:   } else {
2172:     if (diag != 0.0) {
2173:       for (i = 0; i < N; i++) {
2174:         PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2175:         if (a->ilen[rows[i]] > 0) {
2176:           if (rows[i] >= A->cmap->n) {
2177:             a->ilen[rows[i]] = 0;
2178:           } else {
2179:             a->ilen[rows[i]]    = 1;
2180:             aa[a->i[rows[i]]]   = diag;
2181:             a->j[a->i[rows[i]]] = rows[i];
2182:           }
2183:         } else if (rows[i] < A->cmap->n) { /* in case row was completely empty */
2184:           PetscCall(MatSetValues_SeqAIJ(A, 1, &rows[i], 1, &rows[i], &diag, INSERT_VALUES));
2185:         }
2186:       }
2187:     } else {
2188:       for (i = 0; i < N; i++) {
2189:         PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2190:         a->ilen[rows[i]] = 0;
2191:       }
2192:     }
2193:     A->nonzerostate++;
2194:   }
2195:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
2196:   PetscUseTypeMethod(A, assemblyend, MAT_FINAL_ASSEMBLY);
2197:   PetscFunctionReturn(PETSC_SUCCESS);
2198: }

2200: static PetscErrorCode MatZeroRowsColumns_SeqAIJ(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
2201: {
2202:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2203:   PetscInt           i, j, m = A->rmap->n - 1, d = 0;
2204:   PetscBool          missing, *zeroed, vecs = PETSC_FALSE;
2205:   const PetscScalar *xx;
2206:   PetscScalar       *bb, *aa;

2208:   PetscFunctionBegin;
2209:   if (!N) PetscFunctionReturn(PETSC_SUCCESS);
2210:   PetscCall(MatSeqAIJGetArray(A, &aa));
2211:   if (x && b) {
2212:     PetscCall(VecGetArrayRead(x, &xx));
2213:     PetscCall(VecGetArray(b, &bb));
2214:     vecs = PETSC_TRUE;
2215:   }
2216:   PetscCall(PetscCalloc1(A->rmap->n, &zeroed));
2217:   for (i = 0; i < N; i++) {
2218:     PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2219:     PetscCall(PetscArrayzero(PetscSafePointerPlusOffset(aa, a->i[rows[i]]), a->ilen[rows[i]]));

2221:     zeroed[rows[i]] = PETSC_TRUE;
2222:   }
2223:   for (i = 0; i < A->rmap->n; i++) {
2224:     if (!zeroed[i]) {
2225:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
2226:         if (a->j[j] < A->rmap->n && zeroed[a->j[j]]) {
2227:           if (vecs) bb[i] -= aa[j] * xx[a->j[j]];
2228:           aa[j] = 0.0;
2229:         }
2230:       }
2231:     } else if (vecs && i < A->cmap->N) bb[i] = diag * xx[i];
2232:   }
2233:   if (x && b) {
2234:     PetscCall(VecRestoreArrayRead(x, &xx));
2235:     PetscCall(VecRestoreArray(b, &bb));
2236:   }
2237:   PetscCall(PetscFree(zeroed));
2238:   if (diag != 0.0) {
2239:     PetscCall(MatMissingDiagonal_SeqAIJ(A, &missing, &d));
2240:     if (missing) {
2241:       for (i = 0; i < N; i++) {
2242:         if (rows[i] >= A->cmap->N) continue;
2243:         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]);
2244:         PetscCall(MatSetValues_SeqAIJ(A, 1, &rows[i], 1, &rows[i], &diag, INSERT_VALUES));
2245:       }
2246:     } else {
2247:       for (i = 0; i < N; i++) aa[a->diag[rows[i]]] = diag;
2248:     }
2249:   }
2250:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
2251:   PetscUseTypeMethod(A, assemblyend, MAT_FINAL_ASSEMBLY);
2252:   PetscFunctionReturn(PETSC_SUCCESS);
2253: }

2255: PetscErrorCode MatGetRow_SeqAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
2256: {
2257:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2258:   const PetscScalar *aa;

2260:   PetscFunctionBegin;
2261:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2262:   *nz = a->i[row + 1] - a->i[row];
2263:   if (v) *v = PetscSafePointerPlusOffset((PetscScalar *)aa, a->i[row]);
2264:   if (idx) {
2265:     if (*nz && a->j) *idx = a->j + a->i[row];
2266:     else *idx = NULL;
2267:   }
2268:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2269:   PetscFunctionReturn(PETSC_SUCCESS);
2270: }

2272: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
2273: {
2274:   PetscFunctionBegin;
2275:   PetscFunctionReturn(PETSC_SUCCESS);
2276: }

2278: static PetscErrorCode MatNorm_SeqAIJ(Mat A, NormType type, PetscReal *nrm)
2279: {
2280:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
2281:   const MatScalar *v;
2282:   PetscReal        sum = 0.0;
2283:   PetscInt         i, j;

2285:   PetscFunctionBegin;
2286:   PetscCall(MatSeqAIJGetArrayRead(A, &v));
2287:   if (type == NORM_FROBENIUS) {
2288: #if defined(PETSC_USE_REAL___FP16)
2289:     PetscBLASInt one = 1, nz = a->nz;
2290:     PetscCallBLAS("BLASnrm2", *nrm = BLASnrm2_(&nz, v, &one));
2291: #else
2292:     for (i = 0; i < a->nz; i++) {
2293:       sum += PetscRealPart(PetscConj(*v) * (*v));
2294:       v++;
2295:     }
2296:     *nrm = PetscSqrtReal(sum);
2297: #endif
2298:     PetscCall(PetscLogFlops(2.0 * a->nz));
2299:   } else if (type == NORM_1) {
2300:     PetscReal *tmp;
2301:     PetscInt  *jj = a->j;
2302:     PetscCall(PetscCalloc1(A->cmap->n + 1, &tmp));
2303:     *nrm = 0.0;
2304:     for (j = 0; j < a->nz; j++) {
2305:       tmp[*jj++] += PetscAbsScalar(*v);
2306:       v++;
2307:     }
2308:     for (j = 0; j < A->cmap->n; j++) {
2309:       if (tmp[j] > *nrm) *nrm = tmp[j];
2310:     }
2311:     PetscCall(PetscFree(tmp));
2312:     PetscCall(PetscLogFlops(PetscMax(a->nz - 1, 0)));
2313:   } else if (type == NORM_INFINITY) {
2314:     *nrm = 0.0;
2315:     for (j = 0; j < A->rmap->n; j++) {
2316:       const PetscScalar *v2 = PetscSafePointerPlusOffset(v, a->i[j]);
2317:       sum                   = 0.0;
2318:       for (i = 0; i < a->i[j + 1] - a->i[j]; i++) {
2319:         sum += PetscAbsScalar(*v2);
2320:         v2++;
2321:       }
2322:       if (sum > *nrm) *nrm = sum;
2323:     }
2324:     PetscCall(PetscLogFlops(PetscMax(a->nz - 1, 0)));
2325:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for two norm");
2326:   PetscCall(MatSeqAIJRestoreArrayRead(A, &v));
2327:   PetscFunctionReturn(PETSC_SUCCESS);
2328: }

2330: static PetscErrorCode MatIsTranspose_SeqAIJ(Mat A, Mat B, PetscReal tol, PetscBool *f)
2331: {
2332:   Mat_SeqAIJ      *aij = (Mat_SeqAIJ *)A->data, *bij = (Mat_SeqAIJ *)B->data;
2333:   PetscInt        *adx, *bdx, *aii, *bii, *aptr, *bptr;
2334:   const MatScalar *va, *vb;
2335:   PetscInt         ma, na, mb, nb, i;

2337:   PetscFunctionBegin;
2338:   PetscCall(MatGetSize(A, &ma, &na));
2339:   PetscCall(MatGetSize(B, &mb, &nb));
2340:   if (ma != nb || na != mb) {
2341:     *f = PETSC_FALSE;
2342:     PetscFunctionReturn(PETSC_SUCCESS);
2343:   }
2344:   PetscCall(MatSeqAIJGetArrayRead(A, &va));
2345:   PetscCall(MatSeqAIJGetArrayRead(B, &vb));
2346:   aii = aij->i;
2347:   bii = bij->i;
2348:   adx = aij->j;
2349:   bdx = bij->j;
2350:   PetscCall(PetscMalloc1(ma, &aptr));
2351:   PetscCall(PetscMalloc1(mb, &bptr));
2352:   for (i = 0; i < ma; i++) aptr[i] = aii[i];
2353:   for (i = 0; i < mb; i++) bptr[i] = bii[i];

2355:   *f = PETSC_TRUE;
2356:   for (i = 0; i < ma; i++) {
2357:     while (aptr[i] < aii[i + 1]) {
2358:       PetscInt    idc, idr;
2359:       PetscScalar vc, vr;
2360:       /* column/row index/value */
2361:       idc = adx[aptr[i]];
2362:       idr = bdx[bptr[idc]];
2363:       vc  = va[aptr[i]];
2364:       vr  = vb[bptr[idc]];
2365:       if (i != idr || PetscAbsScalar(vc - vr) > tol) {
2366:         *f = PETSC_FALSE;
2367:         goto done;
2368:       } else {
2369:         aptr[i]++;
2370:         if (B || i != idc) bptr[idc]++;
2371:       }
2372:     }
2373:   }
2374: done:
2375:   PetscCall(PetscFree(aptr));
2376:   PetscCall(PetscFree(bptr));
2377:   PetscCall(MatSeqAIJRestoreArrayRead(A, &va));
2378:   PetscCall(MatSeqAIJRestoreArrayRead(B, &vb));
2379:   PetscFunctionReturn(PETSC_SUCCESS);
2380: }

2382: static PetscErrorCode MatIsHermitianTranspose_SeqAIJ(Mat A, Mat B, PetscReal tol, PetscBool *f)
2383: {
2384:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data, *bij = (Mat_SeqAIJ *)B->data;
2385:   PetscInt   *adx, *bdx, *aii, *bii, *aptr, *bptr;
2386:   MatScalar  *va, *vb;
2387:   PetscInt    ma, na, mb, nb, i;

2389:   PetscFunctionBegin;
2390:   PetscCall(MatGetSize(A, &ma, &na));
2391:   PetscCall(MatGetSize(B, &mb, &nb));
2392:   if (ma != nb || na != mb) {
2393:     *f = PETSC_FALSE;
2394:     PetscFunctionReturn(PETSC_SUCCESS);
2395:   }
2396:   aii = aij->i;
2397:   bii = bij->i;
2398:   adx = aij->j;
2399:   bdx = bij->j;
2400:   va  = aij->a;
2401:   vb  = bij->a;
2402:   PetscCall(PetscMalloc1(ma, &aptr));
2403:   PetscCall(PetscMalloc1(mb, &bptr));
2404:   for (i = 0; i < ma; i++) aptr[i] = aii[i];
2405:   for (i = 0; i < mb; i++) bptr[i] = bii[i];

2407:   *f = PETSC_TRUE;
2408:   for (i = 0; i < ma; i++) {
2409:     while (aptr[i] < aii[i + 1]) {
2410:       PetscInt    idc, idr;
2411:       PetscScalar vc, vr;
2412:       /* column/row index/value */
2413:       idc = adx[aptr[i]];
2414:       idr = bdx[bptr[idc]];
2415:       vc  = va[aptr[i]];
2416:       vr  = vb[bptr[idc]];
2417:       if (i != idr || PetscAbsScalar(vc - PetscConj(vr)) > tol) {
2418:         *f = PETSC_FALSE;
2419:         goto done;
2420:       } else {
2421:         aptr[i]++;
2422:         if (B || i != idc) bptr[idc]++;
2423:       }
2424:     }
2425:   }
2426: done:
2427:   PetscCall(PetscFree(aptr));
2428:   PetscCall(PetscFree(bptr));
2429:   PetscFunctionReturn(PETSC_SUCCESS);
2430: }

2432: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A, Vec ll, Vec rr)
2433: {
2434:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2435:   const PetscScalar *l, *r;
2436:   PetscScalar        x;
2437:   MatScalar         *v;
2438:   PetscInt           i, j, m = A->rmap->n, n = A->cmap->n, M, nz = a->nz;
2439:   const PetscInt    *jj;

2441:   PetscFunctionBegin;
2442:   if (ll) {
2443:     /* The local size is used so that VecMPI can be passed to this routine
2444:        by MatDiagonalScale_MPIAIJ */
2445:     PetscCall(VecGetLocalSize(ll, &m));
2446:     PetscCheck(m == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Left scaling vector wrong length");
2447:     PetscCall(VecGetArrayRead(ll, &l));
2448:     PetscCall(MatSeqAIJGetArray(A, &v));
2449:     for (i = 0; i < m; i++) {
2450:       x = l[i];
2451:       M = a->i[i + 1] - a->i[i];
2452:       for (j = 0; j < M; j++) (*v++) *= x;
2453:     }
2454:     PetscCall(VecRestoreArrayRead(ll, &l));
2455:     PetscCall(PetscLogFlops(nz));
2456:     PetscCall(MatSeqAIJRestoreArray(A, &v));
2457:   }
2458:   if (rr) {
2459:     PetscCall(VecGetLocalSize(rr, &n));
2460:     PetscCheck(n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Right scaling vector wrong length");
2461:     PetscCall(VecGetArrayRead(rr, &r));
2462:     PetscCall(MatSeqAIJGetArray(A, &v));
2463:     jj = a->j;
2464:     for (i = 0; i < nz; i++) (*v++) *= r[*jj++];
2465:     PetscCall(MatSeqAIJRestoreArray(A, &v));
2466:     PetscCall(VecRestoreArrayRead(rr, &r));
2467:     PetscCall(PetscLogFlops(nz));
2468:   }
2469:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
2470:   PetscFunctionReturn(PETSC_SUCCESS);
2471: }

2473: PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat A, IS isrow, IS iscol, PetscInt csize, MatReuse scall, Mat *B)
2474: {
2475:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data, *c;
2476:   PetscInt          *smap, i, k, kstart, kend, oldcols = A->cmap->n, *lens;
2477:   PetscInt           row, mat_i, *mat_j, tcol, first, step, *mat_ilen, sum, lensi;
2478:   const PetscInt    *irow, *icol;
2479:   const PetscScalar *aa;
2480:   PetscInt           nrows, ncols;
2481:   PetscInt          *starts, *j_new, *i_new, *aj = a->j, *ai = a->i, ii, *ailen = a->ilen;
2482:   MatScalar         *a_new, *mat_a, *c_a;
2483:   Mat                C;
2484:   PetscBool          stride;

2486:   PetscFunctionBegin;
2487:   PetscCall(ISGetIndices(isrow, &irow));
2488:   PetscCall(ISGetLocalSize(isrow, &nrows));
2489:   PetscCall(ISGetLocalSize(iscol, &ncols));

2491:   PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISSTRIDE, &stride));
2492:   if (stride) {
2493:     PetscCall(ISStrideGetInfo(iscol, &first, &step));
2494:   } else {
2495:     first = 0;
2496:     step  = 0;
2497:   }
2498:   if (stride && step == 1) {
2499:     /* special case of contiguous rows */
2500:     PetscCall(PetscMalloc2(nrows, &lens, nrows, &starts));
2501:     /* loop over new rows determining lens and starting points */
2502:     for (i = 0; i < nrows; i++) {
2503:       kstart    = ai[irow[i]];
2504:       kend      = kstart + ailen[irow[i]];
2505:       starts[i] = kstart;
2506:       for (k = kstart; k < kend; k++) {
2507:         if (aj[k] >= first) {
2508:           starts[i] = k;
2509:           break;
2510:         }
2511:       }
2512:       sum = 0;
2513:       while (k < kend) {
2514:         if (aj[k++] >= first + ncols) break;
2515:         sum++;
2516:       }
2517:       lens[i] = sum;
2518:     }
2519:     /* create submatrix */
2520:     if (scall == MAT_REUSE_MATRIX) {
2521:       PetscInt n_cols, n_rows;
2522:       PetscCall(MatGetSize(*B, &n_rows, &n_cols));
2523:       PetscCheck(n_rows == nrows && n_cols == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reused submatrix wrong size");
2524:       PetscCall(MatZeroEntries(*B));
2525:       C = *B;
2526:     } else {
2527:       PetscInt rbs, cbs;
2528:       PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
2529:       PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE));
2530:       PetscCall(ISGetBlockSize(isrow, &rbs));
2531:       PetscCall(ISGetBlockSize(iscol, &cbs));
2532:       PetscCall(MatSetBlockSizes(C, rbs, cbs));
2533:       PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
2534:       PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(C, 0, lens));
2535:     }
2536:     c = (Mat_SeqAIJ *)C->data;

2538:     /* loop over rows inserting into submatrix */
2539:     PetscCall(MatSeqAIJGetArrayWrite(C, &a_new)); // Not 'a_new = c->a-new', since that raw usage ignores offload state of C
2540:     j_new = c->j;
2541:     i_new = c->i;
2542:     PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2543:     for (i = 0; i < nrows; i++) {
2544:       ii    = starts[i];
2545:       lensi = lens[i];
2546:       if (lensi) {
2547:         for (k = 0; k < lensi; k++) *j_new++ = aj[ii + k] - first;
2548:         PetscCall(PetscArraycpy(a_new, aa + starts[i], lensi));
2549:         a_new += lensi;
2550:       }
2551:       i_new[i + 1] = i_new[i] + lensi;
2552:       c->ilen[i]   = lensi;
2553:     }
2554:     PetscCall(MatSeqAIJRestoreArrayWrite(C, &a_new)); // Set C's offload state properly
2555:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2556:     PetscCall(PetscFree2(lens, starts));
2557:   } else {
2558:     PetscCall(ISGetIndices(iscol, &icol));
2559:     PetscCall(PetscCalloc1(oldcols, &smap));
2560:     PetscCall(PetscMalloc1(1 + nrows, &lens));
2561:     for (i = 0; i < ncols; i++) {
2562:       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);
2563:       smap[icol[i]] = i + 1;
2564:     }

2566:     /* determine lens of each row */
2567:     for (i = 0; i < nrows; i++) {
2568:       kstart  = ai[irow[i]];
2569:       kend    = kstart + a->ilen[irow[i]];
2570:       lens[i] = 0;
2571:       for (k = kstart; k < kend; k++) {
2572:         if (smap[aj[k]]) lens[i]++;
2573:       }
2574:     }
2575:     /* Create and fill new matrix */
2576:     if (scall == MAT_REUSE_MATRIX) {
2577:       PetscBool equal;

2579:       c = (Mat_SeqAIJ *)((*B)->data);
2580:       PetscCheck((*B)->rmap->n == nrows && (*B)->cmap->n == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong size");
2581:       PetscCall(PetscArraycmp(c->ilen, lens, (*B)->rmap->n, &equal));
2582:       PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong number of nonzeros");
2583:       PetscCall(PetscArrayzero(c->ilen, (*B)->rmap->n));
2584:       C = *B;
2585:     } else {
2586:       PetscInt rbs, cbs;
2587:       PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
2588:       PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE));
2589:       PetscCall(ISGetBlockSize(isrow, &rbs));
2590:       PetscCall(ISGetBlockSize(iscol, &cbs));
2591:       if (rbs > 1 || cbs > 1) PetscCall(MatSetBlockSizes(C, rbs, cbs));
2592:       PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
2593:       PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(C, 0, lens));
2594:     }
2595:     PetscCall(MatSeqAIJGetArrayRead(A, &aa));

2597:     c = (Mat_SeqAIJ *)C->data;
2598:     PetscCall(MatSeqAIJGetArrayWrite(C, &c_a)); // Not 'c->a', since that raw usage ignores offload state of C
2599:     for (i = 0; i < nrows; i++) {
2600:       row      = irow[i];
2601:       kstart   = ai[row];
2602:       kend     = kstart + a->ilen[row];
2603:       mat_i    = c->i[i];
2604:       mat_j    = PetscSafePointerPlusOffset(c->j, mat_i);
2605:       mat_a    = PetscSafePointerPlusOffset(c_a, mat_i);
2606:       mat_ilen = c->ilen + i;
2607:       for (k = kstart; k < kend; k++) {
2608:         if ((tcol = smap[a->j[k]])) {
2609:           *mat_j++ = tcol - 1;
2610:           *mat_a++ = aa[k];
2611:           (*mat_ilen)++;
2612:         }
2613:       }
2614:     }
2615:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2616:     /* Free work space */
2617:     PetscCall(ISRestoreIndices(iscol, &icol));
2618:     PetscCall(PetscFree(smap));
2619:     PetscCall(PetscFree(lens));
2620:     /* sort */
2621:     for (i = 0; i < nrows; i++) {
2622:       PetscInt ilen;

2624:       mat_i = c->i[i];
2625:       mat_j = PetscSafePointerPlusOffset(c->j, mat_i);
2626:       mat_a = PetscSafePointerPlusOffset(c_a, mat_i);
2627:       ilen  = c->ilen[i];
2628:       PetscCall(PetscSortIntWithScalarArray(ilen, mat_j, mat_a));
2629:     }
2630:     PetscCall(MatSeqAIJRestoreArrayWrite(C, &c_a));
2631:   }
2632: #if defined(PETSC_HAVE_DEVICE)
2633:   PetscCall(MatBindToCPU(C, A->boundtocpu));
2634: #endif
2635:   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2636:   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));

2638:   PetscCall(ISRestoreIndices(isrow, &irow));
2639:   *B = C;
2640:   PetscFunctionReturn(PETSC_SUCCESS);
2641: }

2643: static PetscErrorCode MatGetMultiProcBlock_SeqAIJ(Mat mat, MPI_Comm subComm, MatReuse scall, Mat *subMat)
2644: {
2645:   Mat B;

2647:   PetscFunctionBegin;
2648:   if (scall == MAT_INITIAL_MATRIX) {
2649:     PetscCall(MatCreate(subComm, &B));
2650:     PetscCall(MatSetSizes(B, mat->rmap->n, mat->cmap->n, mat->rmap->n, mat->cmap->n));
2651:     PetscCall(MatSetBlockSizesFromMats(B, mat, mat));
2652:     PetscCall(MatSetType(B, MATSEQAIJ));
2653:     PetscCall(MatDuplicateNoCreate_SeqAIJ(B, mat, MAT_COPY_VALUES, PETSC_TRUE));
2654:     *subMat = B;
2655:   } else {
2656:     PetscCall(MatCopy_SeqAIJ(mat, *subMat, SAME_NONZERO_PATTERN));
2657:   }
2658:   PetscFunctionReturn(PETSC_SUCCESS);
2659: }

2661: static PetscErrorCode MatILUFactor_SeqAIJ(Mat inA, IS row, IS col, const MatFactorInfo *info)
2662: {
2663:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)inA->data;
2664:   Mat         outA;
2665:   PetscBool   row_identity, col_identity;

2667:   PetscFunctionBegin;
2668:   PetscCheck(info->levels == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only levels=0 supported for in-place ilu");

2670:   PetscCall(ISIdentity(row, &row_identity));
2671:   PetscCall(ISIdentity(col, &col_identity));

2673:   outA             = inA;
2674:   outA->factortype = MAT_FACTOR_LU;
2675:   PetscCall(PetscFree(inA->solvertype));
2676:   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &inA->solvertype));

2678:   PetscCall(PetscObjectReference((PetscObject)row));
2679:   PetscCall(ISDestroy(&a->row));

2681:   a->row = row;

2683:   PetscCall(PetscObjectReference((PetscObject)col));
2684:   PetscCall(ISDestroy(&a->col));

2686:   a->col = col;

2688:   /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
2689:   PetscCall(ISDestroy(&a->icol));
2690:   PetscCall(ISInvertPermutation(col, PETSC_DECIDE, &a->icol));

2692:   if (!a->solve_work) { /* this matrix may have been factored before */
2693:     PetscCall(PetscMalloc1(inA->rmap->n + 1, &a->solve_work));
2694:   }

2696:   PetscCall(MatMarkDiagonal_SeqAIJ(inA));
2697:   if (row_identity && col_identity) {
2698:     PetscCall(MatLUFactorNumeric_SeqAIJ_inplace(outA, inA, info));
2699:   } else {
2700:     PetscCall(MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(outA, inA, info));
2701:   }
2702:   PetscFunctionReturn(PETSC_SUCCESS);
2703: }

2705: PetscErrorCode MatScale_SeqAIJ(Mat inA, PetscScalar alpha)
2706: {
2707:   Mat_SeqAIJ  *a = (Mat_SeqAIJ *)inA->data;
2708:   PetscScalar *v;
2709:   PetscBLASInt one = 1, bnz;

2711:   PetscFunctionBegin;
2712:   PetscCall(MatSeqAIJGetArray(inA, &v));
2713:   PetscCall(PetscBLASIntCast(a->nz, &bnz));
2714:   PetscCallBLAS("BLASscal", BLASscal_(&bnz, &alpha, v, &one));
2715:   PetscCall(PetscLogFlops(a->nz));
2716:   PetscCall(MatSeqAIJRestoreArray(inA, &v));
2717:   PetscCall(MatSeqAIJInvalidateDiagonal(inA));
2718:   PetscFunctionReturn(PETSC_SUCCESS);
2719: }

2721: PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *submatj)
2722: {
2723:   PetscInt i;

2725:   PetscFunctionBegin;
2726:   if (!submatj->id) { /* delete data that are linked only to submats[id=0] */
2727:     PetscCall(PetscFree4(submatj->sbuf1, submatj->ptr, submatj->tmp, submatj->ctr));

2729:     for (i = 0; i < submatj->nrqr; ++i) PetscCall(PetscFree(submatj->sbuf2[i]));
2730:     PetscCall(PetscFree3(submatj->sbuf2, submatj->req_size, submatj->req_source1));

2732:     if (submatj->rbuf1) {
2733:       PetscCall(PetscFree(submatj->rbuf1[0]));
2734:       PetscCall(PetscFree(submatj->rbuf1));
2735:     }

2737:     for (i = 0; i < submatj->nrqs; ++i) PetscCall(PetscFree(submatj->rbuf3[i]));
2738:     PetscCall(PetscFree3(submatj->req_source2, submatj->rbuf2, submatj->rbuf3));
2739:     PetscCall(PetscFree(submatj->pa));
2740:   }

2742: #if defined(PETSC_USE_CTABLE)
2743:   PetscCall(PetscHMapIDestroy(&submatj->rmap));
2744:   if (submatj->cmap_loc) PetscCall(PetscFree(submatj->cmap_loc));
2745:   PetscCall(PetscFree(submatj->rmap_loc));
2746: #else
2747:   PetscCall(PetscFree(submatj->rmap));
2748: #endif

2750:   if (!submatj->allcolumns) {
2751: #if defined(PETSC_USE_CTABLE)
2752:     PetscCall(PetscHMapIDestroy(&submatj->cmap));
2753: #else
2754:     PetscCall(PetscFree(submatj->cmap));
2755: #endif
2756:   }
2757:   PetscCall(PetscFree(submatj->row2proc));

2759:   PetscCall(PetscFree(submatj));
2760:   PetscFunctionReturn(PETSC_SUCCESS);
2761: }

2763: PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat C)
2764: {
2765:   Mat_SeqAIJ  *c       = (Mat_SeqAIJ *)C->data;
2766:   Mat_SubSppt *submatj = c->submatis1;

2768:   PetscFunctionBegin;
2769:   PetscCall((*submatj->destroy)(C));
2770:   PetscCall(MatDestroySubMatrix_Private(submatj));
2771:   PetscFunctionReturn(PETSC_SUCCESS);
2772: }

2774: /* Note this has code duplication with MatDestroySubMatrices_SeqBAIJ() */
2775: static PetscErrorCode MatDestroySubMatrices_SeqAIJ(PetscInt n, Mat *mat[])
2776: {
2777:   PetscInt     i;
2778:   Mat          C;
2779:   Mat_SeqAIJ  *c;
2780:   Mat_SubSppt *submatj;

2782:   PetscFunctionBegin;
2783:   for (i = 0; i < n; i++) {
2784:     C       = (*mat)[i];
2785:     c       = (Mat_SeqAIJ *)C->data;
2786:     submatj = c->submatis1;
2787:     if (submatj) {
2788:       if (--((PetscObject)C)->refct <= 0) {
2789:         PetscCall(PetscFree(C->factorprefix));
2790:         PetscCall((*submatj->destroy)(C));
2791:         PetscCall(MatDestroySubMatrix_Private(submatj));
2792:         PetscCall(PetscFree(C->defaultvectype));
2793:         PetscCall(PetscFree(C->defaultrandtype));
2794:         PetscCall(PetscLayoutDestroy(&C->rmap));
2795:         PetscCall(PetscLayoutDestroy(&C->cmap));
2796:         PetscCall(PetscHeaderDestroy(&C));
2797:       }
2798:     } else {
2799:       PetscCall(MatDestroy(&C));
2800:     }
2801:   }

2803:   /* Destroy Dummy submatrices created for reuse */
2804:   PetscCall(MatDestroySubMatrices_Dummy(n, mat));

2806:   PetscCall(PetscFree(*mat));
2807:   PetscFunctionReturn(PETSC_SUCCESS);
2808: }

2810: static PetscErrorCode MatCreateSubMatrices_SeqAIJ(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *B[])
2811: {
2812:   PetscInt i;

2814:   PetscFunctionBegin;
2815:   if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscCalloc1(n + 1, B));

2817:   for (i = 0; i < n; i++) PetscCall(MatCreateSubMatrix_SeqAIJ(A, irow[i], icol[i], PETSC_DECIDE, scall, &(*B)[i]));
2818:   PetscFunctionReturn(PETSC_SUCCESS);
2819: }

2821: static PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A, PetscInt is_max, IS is[], PetscInt ov)
2822: {
2823:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
2824:   PetscInt        row, i, j, k, l, ll, m, n, *nidx, isz, val;
2825:   const PetscInt *idx;
2826:   PetscInt        start, end, *ai, *aj, bs = (A->rmap->bs > 0 && A->rmap->bs == A->cmap->bs) ? A->rmap->bs : 1;
2827:   PetscBT         table;

2829:   PetscFunctionBegin;
2830:   m  = A->rmap->n / bs;
2831:   ai = a->i;
2832:   aj = a->j;

2834:   PetscCheck(ov >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "illegal negative overlap value used");

2836:   PetscCall(PetscMalloc1(m + 1, &nidx));
2837:   PetscCall(PetscBTCreate(m, &table));

2839:   for (i = 0; i < is_max; i++) {
2840:     /* Initialize the two local arrays */
2841:     isz = 0;
2842:     PetscCall(PetscBTMemzero(m, table));

2844:     /* Extract the indices, assume there can be duplicate entries */
2845:     PetscCall(ISGetIndices(is[i], &idx));
2846:     PetscCall(ISGetLocalSize(is[i], &n));

2848:     if (bs > 1) {
2849:       /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
2850:       for (j = 0; j < n; ++j) {
2851:         if (!PetscBTLookupSet(table, idx[j] / bs)) nidx[isz++] = idx[j] / bs;
2852:       }
2853:       PetscCall(ISRestoreIndices(is[i], &idx));
2854:       PetscCall(ISDestroy(&is[i]));

2856:       k = 0;
2857:       for (j = 0; j < ov; j++) { /* for each overlap */
2858:         n = isz;
2859:         for (; k < n; k++) { /* do only those rows in nidx[k], which are not done yet */
2860:           for (ll = 0; ll < bs; ll++) {
2861:             row   = bs * nidx[k] + ll;
2862:             start = ai[row];
2863:             end   = ai[row + 1];
2864:             for (l = start; l < end; l++) {
2865:               val = aj[l] / bs;
2866:               if (!PetscBTLookupSet(table, val)) nidx[isz++] = val;
2867:             }
2868:           }
2869:         }
2870:       }
2871:       PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, isz, nidx, PETSC_COPY_VALUES, is + i));
2872:     } else {
2873:       /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
2874:       for (j = 0; j < n; ++j) {
2875:         if (!PetscBTLookupSet(table, idx[j])) nidx[isz++] = idx[j];
2876:       }
2877:       PetscCall(ISRestoreIndices(is[i], &idx));
2878:       PetscCall(ISDestroy(&is[i]));

2880:       k = 0;
2881:       for (j = 0; j < ov; j++) { /* for each overlap */
2882:         n = isz;
2883:         for (; k < n; k++) { /* do only those rows in nidx[k], which are not done yet */
2884:           row   = nidx[k];
2885:           start = ai[row];
2886:           end   = ai[row + 1];
2887:           for (l = start; l < end; l++) {
2888:             val = aj[l];
2889:             if (!PetscBTLookupSet(table, val)) nidx[isz++] = val;
2890:           }
2891:         }
2892:       }
2893:       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, PETSC_COPY_VALUES, is + i));
2894:     }
2895:   }
2896:   PetscCall(PetscBTDestroy(&table));
2897:   PetscCall(PetscFree(nidx));
2898:   PetscFunctionReturn(PETSC_SUCCESS);
2899: }

2901: static PetscErrorCode MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
2902: {
2903:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
2904:   PetscInt        i, nz = 0, m = A->rmap->n, n = A->cmap->n;
2905:   const PetscInt *row, *col;
2906:   PetscInt       *cnew, j, *lens;
2907:   IS              icolp, irowp;
2908:   PetscInt       *cwork = NULL;
2909:   PetscScalar    *vwork = NULL;

2911:   PetscFunctionBegin;
2912:   PetscCall(ISInvertPermutation(rowp, PETSC_DECIDE, &irowp));
2913:   PetscCall(ISGetIndices(irowp, &row));
2914:   PetscCall(ISInvertPermutation(colp, PETSC_DECIDE, &icolp));
2915:   PetscCall(ISGetIndices(icolp, &col));

2917:   /* determine lengths of permuted rows */
2918:   PetscCall(PetscMalloc1(m + 1, &lens));
2919:   for (i = 0; i < m; i++) lens[row[i]] = a->i[i + 1] - a->i[i];
2920:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
2921:   PetscCall(MatSetSizes(*B, m, n, m, n));
2922:   PetscCall(MatSetBlockSizesFromMats(*B, A, A));
2923:   PetscCall(MatSetType(*B, ((PetscObject)A)->type_name));
2924:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*B, 0, lens));
2925:   PetscCall(PetscFree(lens));

2927:   PetscCall(PetscMalloc1(n, &cnew));
2928:   for (i = 0; i < m; i++) {
2929:     PetscCall(MatGetRow_SeqAIJ(A, i, &nz, &cwork, &vwork));
2930:     for (j = 0; j < nz; j++) cnew[j] = col[cwork[j]];
2931:     PetscCall(MatSetValues_SeqAIJ(*B, 1, &row[i], nz, cnew, vwork, INSERT_VALUES));
2932:     PetscCall(MatRestoreRow_SeqAIJ(A, i, &nz, &cwork, &vwork));
2933:   }
2934:   PetscCall(PetscFree(cnew));

2936:   (*B)->assembled = PETSC_FALSE;

2938: #if defined(PETSC_HAVE_DEVICE)
2939:   PetscCall(MatBindToCPU(*B, A->boundtocpu));
2940: #endif
2941:   PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY));
2942:   PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY));
2943:   PetscCall(ISRestoreIndices(irowp, &row));
2944:   PetscCall(ISRestoreIndices(icolp, &col));
2945:   PetscCall(ISDestroy(&irowp));
2946:   PetscCall(ISDestroy(&icolp));
2947:   if (rowp == colp) PetscCall(MatPropagateSymmetryOptions(A, *B));
2948:   PetscFunctionReturn(PETSC_SUCCESS);
2949: }

2951: PetscErrorCode MatCopy_SeqAIJ(Mat A, Mat B, MatStructure str)
2952: {
2953:   PetscFunctionBegin;
2954:   /* If the two matrices have the same copy implementation, use fast copy. */
2955:   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
2956:     Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2957:     Mat_SeqAIJ        *b = (Mat_SeqAIJ *)B->data;
2958:     const PetscScalar *aa;

2960:     PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2961:     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]);
2962:     PetscCall(PetscArraycpy(b->a, aa, a->i[A->rmap->n]));
2963:     PetscCall(PetscObjectStateIncrease((PetscObject)B));
2964:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2965:   } else {
2966:     PetscCall(MatCopy_Basic(A, B, str));
2967:   }
2968:   PetscFunctionReturn(PETSC_SUCCESS);
2969: }

2971: PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat A, PetscScalar *array[])
2972: {
2973:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

2975:   PetscFunctionBegin;
2976:   *array = a->a;
2977:   PetscFunctionReturn(PETSC_SUCCESS);
2978: }

2980: PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat A, PetscScalar *array[])
2981: {
2982:   PetscFunctionBegin;
2983:   *array = NULL;
2984:   PetscFunctionReturn(PETSC_SUCCESS);
2985: }

2987: /*
2988:    Computes the number of nonzeros per row needed for preallocation when X and Y
2989:    have different nonzero structure.
2990: */
2991: PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt m, const PetscInt *xi, const PetscInt *xj, const PetscInt *yi, const PetscInt *yj, PetscInt *nnz)
2992: {
2993:   PetscInt i, j, k, nzx, nzy;

2995:   PetscFunctionBegin;
2996:   /* Set the number of nonzeros in the new matrix */
2997:   for (i = 0; i < m; i++) {
2998:     const PetscInt *xjj = PetscSafePointerPlusOffset(xj, xi[i]), *yjj = PetscSafePointerPlusOffset(yj, yi[i]);
2999:     nzx    = xi[i + 1] - xi[i];
3000:     nzy    = yi[i + 1] - yi[i];
3001:     nnz[i] = 0;
3002:     for (j = 0, k = 0; j < nzx; j++) {                  /* Point in X */
3003:       for (; k < nzy && yjj[k] < xjj[j]; k++) nnz[i]++; /* Catch up to X */
3004:       if (k < nzy && yjj[k] == xjj[j]) k++;             /* Skip duplicate */
3005:       nnz[i]++;
3006:     }
3007:     for (; k < nzy; k++) nnz[i]++;
3008:   }
3009:   PetscFunctionReturn(PETSC_SUCCESS);
3010: }

3012: PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat Y, Mat X, PetscInt *nnz)
3013: {
3014:   PetscInt    m = Y->rmap->N;
3015:   Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data;
3016:   Mat_SeqAIJ *y = (Mat_SeqAIJ *)Y->data;

3018:   PetscFunctionBegin;
3019:   /* Set the number of nonzeros in the new matrix */
3020:   PetscCall(MatAXPYGetPreallocation_SeqX_private(m, x->i, x->j, y->i, y->j, nnz));
3021:   PetscFunctionReturn(PETSC_SUCCESS);
3022: }

3024: PetscErrorCode MatAXPY_SeqAIJ(Mat Y, PetscScalar a, Mat X, MatStructure str)
3025: {
3026:   Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data, *y = (Mat_SeqAIJ *)Y->data;

3028:   PetscFunctionBegin;
3029:   if (str == UNKNOWN_NONZERO_PATTERN || (PetscDefined(USE_DEBUG) && str == SAME_NONZERO_PATTERN)) {
3030:     PetscBool e = x->nz == y->nz ? PETSC_TRUE : PETSC_FALSE;
3031:     if (e) {
3032:       PetscCall(PetscArraycmp(x->i, y->i, Y->rmap->n + 1, &e));
3033:       if (e) {
3034:         PetscCall(PetscArraycmp(x->j, y->j, y->nz, &e));
3035:         if (e) str = SAME_NONZERO_PATTERN;
3036:       }
3037:     }
3038:     if (!e) PetscCheck(str != SAME_NONZERO_PATTERN, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatStructure is not SAME_NONZERO_PATTERN");
3039:   }
3040:   if (str == SAME_NONZERO_PATTERN) {
3041:     const PetscScalar *xa;
3042:     PetscScalar       *ya, alpha = a;
3043:     PetscBLASInt       one = 1, bnz;

3045:     PetscCall(PetscBLASIntCast(x->nz, &bnz));
3046:     PetscCall(MatSeqAIJGetArray(Y, &ya));
3047:     PetscCall(MatSeqAIJGetArrayRead(X, &xa));
3048:     PetscCallBLAS("BLASaxpy", BLASaxpy_(&bnz, &alpha, xa, &one, ya, &one));
3049:     PetscCall(MatSeqAIJRestoreArrayRead(X, &xa));
3050:     PetscCall(MatSeqAIJRestoreArray(Y, &ya));
3051:     PetscCall(PetscLogFlops(2.0 * bnz));
3052:     PetscCall(MatSeqAIJInvalidateDiagonal(Y));
3053:     PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3054:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
3055:     PetscCall(MatAXPY_Basic(Y, a, X, str));
3056:   } else {
3057:     Mat       B;
3058:     PetscInt *nnz;
3059:     PetscCall(PetscMalloc1(Y->rmap->N, &nnz));
3060:     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &B));
3061:     PetscCall(PetscObjectSetName((PetscObject)B, ((PetscObject)Y)->name));
3062:     PetscCall(MatSetLayouts(B, Y->rmap, Y->cmap));
3063:     PetscCall(MatSetType(B, ((PetscObject)Y)->type_name));
3064:     PetscCall(MatAXPYGetPreallocation_SeqAIJ(Y, X, nnz));
3065:     PetscCall(MatSeqAIJSetPreallocation(B, 0, nnz));
3066:     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
3067:     PetscCall(MatHeaderMerge(Y, &B));
3068:     PetscCall(MatSeqAIJCheckInode(Y));
3069:     PetscCall(PetscFree(nnz));
3070:   }
3071:   PetscFunctionReturn(PETSC_SUCCESS);
3072: }

3074: PETSC_INTERN PetscErrorCode MatConjugate_SeqAIJ(Mat mat)
3075: {
3076: #if defined(PETSC_USE_COMPLEX)
3077:   Mat_SeqAIJ  *aij = (Mat_SeqAIJ *)mat->data;
3078:   PetscInt     i, nz;
3079:   PetscScalar *a;

3081:   PetscFunctionBegin;
3082:   nz = aij->nz;
3083:   PetscCall(MatSeqAIJGetArray(mat, &a));
3084:   for (i = 0; i < nz; i++) a[i] = PetscConj(a[i]);
3085:   PetscCall(MatSeqAIJRestoreArray(mat, &a));
3086: #else
3087:   PetscFunctionBegin;
3088: #endif
3089:   PetscFunctionReturn(PETSC_SUCCESS);
3090: }

3092: static PetscErrorCode MatGetRowMaxAbs_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3093: {
3094:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3095:   PetscInt         i, j, m = A->rmap->n, *ai, *aj, ncols, n;
3096:   PetscReal        atmp;
3097:   PetscScalar     *x;
3098:   const MatScalar *aa, *av;

3100:   PetscFunctionBegin;
3101:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3102:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3103:   aa = av;
3104:   ai = a->i;
3105:   aj = a->j;

3107:   PetscCall(VecSet(v, 0.0));
3108:   PetscCall(VecGetArrayWrite(v, &x));
3109:   PetscCall(VecGetLocalSize(v, &n));
3110:   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3111:   for (i = 0; i < m; i++) {
3112:     ncols = ai[1] - ai[0];
3113:     ai++;
3114:     for (j = 0; j < ncols; j++) {
3115:       atmp = PetscAbsScalar(*aa);
3116:       if (PetscAbsScalar(x[i]) < atmp) {
3117:         x[i] = atmp;
3118:         if (idx) idx[i] = *aj;
3119:       }
3120:       aa++;
3121:       aj++;
3122:     }
3123:   }
3124:   PetscCall(VecRestoreArrayWrite(v, &x));
3125:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3126:   PetscFunctionReturn(PETSC_SUCCESS);
3127: }

3129: static PetscErrorCode MatGetRowSumAbs_SeqAIJ(Mat A, Vec v)
3130: {
3131:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3132:   PetscInt         i, j, m = A->rmap->n, *ai, ncols, n;
3133:   PetscScalar     *x;
3134:   const MatScalar *aa, *av;

3136:   PetscFunctionBegin;
3137:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3138:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3139:   aa = av;
3140:   ai = a->i;

3142:   PetscCall(VecSet(v, 0.0));
3143:   PetscCall(VecGetArrayWrite(v, &x));
3144:   PetscCall(VecGetLocalSize(v, &n));
3145:   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3146:   for (i = 0; i < m; i++) {
3147:     ncols = ai[1] - ai[0];
3148:     ai++;
3149:     for (j = 0; j < ncols; j++) {
3150:       x[i] += PetscAbsScalar(*aa);
3151:       aa++;
3152:     }
3153:   }
3154:   PetscCall(VecRestoreArrayWrite(v, &x));
3155:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3156:   PetscFunctionReturn(PETSC_SUCCESS);
3157: }

3159: static PetscErrorCode MatGetRowMax_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3160: {
3161:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3162:   PetscInt         i, j, m = A->rmap->n, *ai, *aj, ncols, n;
3163:   PetscScalar     *x;
3164:   const MatScalar *aa, *av;

3166:   PetscFunctionBegin;
3167:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3168:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3169:   aa = av;
3170:   ai = a->i;
3171:   aj = a->j;

3173:   PetscCall(VecSet(v, 0.0));
3174:   PetscCall(VecGetArrayWrite(v, &x));
3175:   PetscCall(VecGetLocalSize(v, &n));
3176:   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3177:   for (i = 0; i < m; i++) {
3178:     ncols = ai[1] - ai[0];
3179:     ai++;
3180:     if (ncols == A->cmap->n) { /* row is dense */
3181:       x[i] = *aa;
3182:       if (idx) idx[i] = 0;
3183:     } else { /* row is sparse so already KNOW maximum is 0.0 or higher */
3184:       x[i] = 0.0;
3185:       if (idx) {
3186:         for (j = 0; j < ncols; j++) { /* find first implicit 0.0 in the row */
3187:           if (aj[j] > j) {
3188:             idx[i] = j;
3189:             break;
3190:           }
3191:         }
3192:         /* in case first implicit 0.0 in the row occurs at ncols-th column */
3193:         if (j == ncols && j < A->cmap->n) idx[i] = j;
3194:       }
3195:     }
3196:     for (j = 0; j < ncols; j++) {
3197:       if (PetscRealPart(x[i]) < PetscRealPart(*aa)) {
3198:         x[i] = *aa;
3199:         if (idx) idx[i] = *aj;
3200:       }
3201:       aa++;
3202:       aj++;
3203:     }
3204:   }
3205:   PetscCall(VecRestoreArrayWrite(v, &x));
3206:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3207:   PetscFunctionReturn(PETSC_SUCCESS);
3208: }

3210: static PetscErrorCode MatGetRowMinAbs_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3211: {
3212:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3213:   PetscInt         i, j, m = A->rmap->n, *ai, *aj, ncols, n;
3214:   PetscScalar     *x;
3215:   const MatScalar *aa, *av;

3217:   PetscFunctionBegin;
3218:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3219:   aa = av;
3220:   ai = a->i;
3221:   aj = a->j;

3223:   PetscCall(VecSet(v, 0.0));
3224:   PetscCall(VecGetArrayWrite(v, &x));
3225:   PetscCall(VecGetLocalSize(v, &n));
3226:   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector, %" PetscInt_FMT " vs. %" PetscInt_FMT " rows", m, n);
3227:   for (i = 0; i < m; i++) {
3228:     ncols = ai[1] - ai[0];
3229:     ai++;
3230:     if (ncols == A->cmap->n) { /* row is dense */
3231:       x[i] = *aa;
3232:       if (idx) idx[i] = 0;
3233:     } else { /* row is sparse so already KNOW minimum is 0.0 or higher */
3234:       x[i] = 0.0;
3235:       if (idx) { /* find first implicit 0.0 in the row */
3236:         for (j = 0; j < ncols; j++) {
3237:           if (aj[j] > j) {
3238:             idx[i] = j;
3239:             break;
3240:           }
3241:         }
3242:         /* in case first implicit 0.0 in the row occurs at ncols-th column */
3243:         if (j == ncols && j < A->cmap->n) idx[i] = j;
3244:       }
3245:     }
3246:     for (j = 0; j < ncols; j++) {
3247:       if (PetscAbsScalar(x[i]) > PetscAbsScalar(*aa)) {
3248:         x[i] = *aa;
3249:         if (idx) idx[i] = *aj;
3250:       }
3251:       aa++;
3252:       aj++;
3253:     }
3254:   }
3255:   PetscCall(VecRestoreArrayWrite(v, &x));
3256:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3257:   PetscFunctionReturn(PETSC_SUCCESS);
3258: }

3260: static PetscErrorCode MatGetRowMin_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3261: {
3262:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3263:   PetscInt         i, j, m = A->rmap->n, ncols, n;
3264:   const PetscInt  *ai, *aj;
3265:   PetscScalar     *x;
3266:   const MatScalar *aa, *av;

3268:   PetscFunctionBegin;
3269:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3270:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3271:   aa = av;
3272:   ai = a->i;
3273:   aj = a->j;

3275:   PetscCall(VecSet(v, 0.0));
3276:   PetscCall(VecGetArrayWrite(v, &x));
3277:   PetscCall(VecGetLocalSize(v, &n));
3278:   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3279:   for (i = 0; i < m; i++) {
3280:     ncols = ai[1] - ai[0];
3281:     ai++;
3282:     if (ncols == A->cmap->n) { /* row is dense */
3283:       x[i] = *aa;
3284:       if (idx) idx[i] = 0;
3285:     } else { /* row is sparse so already KNOW minimum is 0.0 or lower */
3286:       x[i] = 0.0;
3287:       if (idx) { /* find first implicit 0.0 in the row */
3288:         for (j = 0; j < ncols; j++) {
3289:           if (aj[j] > j) {
3290:             idx[i] = j;
3291:             break;
3292:           }
3293:         }
3294:         /* in case first implicit 0.0 in the row occurs at ncols-th column */
3295:         if (j == ncols && j < A->cmap->n) idx[i] = j;
3296:       }
3297:     }
3298:     for (j = 0; j < ncols; j++) {
3299:       if (PetscRealPart(x[i]) > PetscRealPart(*aa)) {
3300:         x[i] = *aa;
3301:         if (idx) idx[i] = *aj;
3302:       }
3303:       aa++;
3304:       aj++;
3305:     }
3306:   }
3307:   PetscCall(VecRestoreArrayWrite(v, &x));
3308:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3309:   PetscFunctionReturn(PETSC_SUCCESS);
3310: }

3312: static PetscErrorCode MatInvertBlockDiagonal_SeqAIJ(Mat A, const PetscScalar **values)
3313: {
3314:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
3315:   PetscInt        i, bs = PetscAbs(A->rmap->bs), mbs = A->rmap->n / bs, ipvt[5], bs2 = bs * bs, *v_pivots, ij[7], *IJ, j;
3316:   MatScalar      *diag, work[25], *v_work;
3317:   const PetscReal shift = 0.0;
3318:   PetscBool       allowzeropivot, zeropivotdetected = PETSC_FALSE;

3320:   PetscFunctionBegin;
3321:   allowzeropivot = PetscNot(A->erroriffailure);
3322:   if (a->ibdiagvalid) {
3323:     if (values) *values = a->ibdiag;
3324:     PetscFunctionReturn(PETSC_SUCCESS);
3325:   }
3326:   PetscCall(MatMarkDiagonal_SeqAIJ(A));
3327:   if (!a->ibdiag) { PetscCall(PetscMalloc1(bs2 * mbs, &a->ibdiag)); }
3328:   diag = a->ibdiag;
3329:   if (values) *values = a->ibdiag;
3330:   /* factor and invert each block */
3331:   switch (bs) {
3332:   case 1:
3333:     for (i = 0; i < mbs; i++) {
3334:       PetscCall(MatGetValues(A, 1, &i, 1, &i, diag + i));
3335:       if (PetscAbsScalar(diag[i] + shift) < PETSC_MACHINE_EPSILON) {
3336:         if (allowzeropivot) {
3337:           A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3338:           A->factorerror_zeropivot_value = PetscAbsScalar(diag[i]);
3339:           A->factorerror_zeropivot_row   = i;
3340:           PetscCall(PetscInfo(A, "Zero pivot, row %" PetscInt_FMT " pivot %g tolerance %g\n", i, (double)PetscAbsScalar(diag[i]), (double)PETSC_MACHINE_EPSILON));
3341:         } 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);
3342:       }
3343:       diag[i] = (PetscScalar)1.0 / (diag[i] + shift);
3344:     }
3345:     break;
3346:   case 2:
3347:     for (i = 0; i < mbs; i++) {
3348:       ij[0] = 2 * i;
3349:       ij[1] = 2 * i + 1;
3350:       PetscCall(MatGetValues(A, 2, ij, 2, ij, diag));
3351:       PetscCall(PetscKernel_A_gets_inverse_A_2(diag, shift, allowzeropivot, &zeropivotdetected));
3352:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3353:       PetscCall(PetscKernel_A_gets_transpose_A_2(diag));
3354:       diag += 4;
3355:     }
3356:     break;
3357:   case 3:
3358:     for (i = 0; i < mbs; i++) {
3359:       ij[0] = 3 * i;
3360:       ij[1] = 3 * i + 1;
3361:       ij[2] = 3 * i + 2;
3362:       PetscCall(MatGetValues(A, 3, ij, 3, ij, diag));
3363:       PetscCall(PetscKernel_A_gets_inverse_A_3(diag, shift, allowzeropivot, &zeropivotdetected));
3364:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3365:       PetscCall(PetscKernel_A_gets_transpose_A_3(diag));
3366:       diag += 9;
3367:     }
3368:     break;
3369:   case 4:
3370:     for (i = 0; i < mbs; i++) {
3371:       ij[0] = 4 * i;
3372:       ij[1] = 4 * i + 1;
3373:       ij[2] = 4 * i + 2;
3374:       ij[3] = 4 * i + 3;
3375:       PetscCall(MatGetValues(A, 4, ij, 4, ij, diag));
3376:       PetscCall(PetscKernel_A_gets_inverse_A_4(diag, shift, allowzeropivot, &zeropivotdetected));
3377:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3378:       PetscCall(PetscKernel_A_gets_transpose_A_4(diag));
3379:       diag += 16;
3380:     }
3381:     break;
3382:   case 5:
3383:     for (i = 0; i < mbs; i++) {
3384:       ij[0] = 5 * i;
3385:       ij[1] = 5 * i + 1;
3386:       ij[2] = 5 * i + 2;
3387:       ij[3] = 5 * i + 3;
3388:       ij[4] = 5 * i + 4;
3389:       PetscCall(MatGetValues(A, 5, ij, 5, ij, diag));
3390:       PetscCall(PetscKernel_A_gets_inverse_A_5(diag, ipvt, work, shift, allowzeropivot, &zeropivotdetected));
3391:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3392:       PetscCall(PetscKernel_A_gets_transpose_A_5(diag));
3393:       diag += 25;
3394:     }
3395:     break;
3396:   case 6:
3397:     for (i = 0; i < mbs; i++) {
3398:       ij[0] = 6 * i;
3399:       ij[1] = 6 * i + 1;
3400:       ij[2] = 6 * i + 2;
3401:       ij[3] = 6 * i + 3;
3402:       ij[4] = 6 * i + 4;
3403:       ij[5] = 6 * i + 5;
3404:       PetscCall(MatGetValues(A, 6, ij, 6, ij, diag));
3405:       PetscCall(PetscKernel_A_gets_inverse_A_6(diag, shift, allowzeropivot, &zeropivotdetected));
3406:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3407:       PetscCall(PetscKernel_A_gets_transpose_A_6(diag));
3408:       diag += 36;
3409:     }
3410:     break;
3411:   case 7:
3412:     for (i = 0; i < mbs; i++) {
3413:       ij[0] = 7 * i;
3414:       ij[1] = 7 * i + 1;
3415:       ij[2] = 7 * i + 2;
3416:       ij[3] = 7 * i + 3;
3417:       ij[4] = 7 * i + 4;
3418:       ij[5] = 7 * i + 5;
3419:       ij[6] = 7 * i + 6;
3420:       PetscCall(MatGetValues(A, 7, ij, 7, ij, diag));
3421:       PetscCall(PetscKernel_A_gets_inverse_A_7(diag, shift, allowzeropivot, &zeropivotdetected));
3422:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3423:       PetscCall(PetscKernel_A_gets_transpose_A_7(diag));
3424:       diag += 49;
3425:     }
3426:     break;
3427:   default:
3428:     PetscCall(PetscMalloc3(bs, &v_work, bs, &v_pivots, bs, &IJ));
3429:     for (i = 0; i < mbs; i++) {
3430:       for (j = 0; j < bs; j++) IJ[j] = bs * i + j;
3431:       PetscCall(MatGetValues(A, bs, IJ, bs, IJ, diag));
3432:       PetscCall(PetscKernel_A_gets_inverse_A(bs, diag, v_pivots, v_work, allowzeropivot, &zeropivotdetected));
3433:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3434:       PetscCall(PetscKernel_A_gets_transpose_A_N(diag, bs));
3435:       diag += bs2;
3436:     }
3437:     PetscCall(PetscFree3(v_work, v_pivots, IJ));
3438:   }
3439:   a->ibdiagvalid = PETSC_TRUE;
3440:   PetscFunctionReturn(PETSC_SUCCESS);
3441: }

3443: static PetscErrorCode MatSetRandom_SeqAIJ(Mat x, PetscRandom rctx)
3444: {
3445:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)x->data;
3446:   PetscScalar a, *aa;
3447:   PetscInt    m, n, i, j, col;

3449:   PetscFunctionBegin;
3450:   if (!x->assembled) {
3451:     PetscCall(MatGetSize(x, &m, &n));
3452:     for (i = 0; i < m; i++) {
3453:       for (j = 0; j < aij->imax[i]; j++) {
3454:         PetscCall(PetscRandomGetValue(rctx, &a));
3455:         col = (PetscInt)(n * PetscRealPart(a));
3456:         PetscCall(MatSetValues(x, 1, &i, 1, &col, &a, ADD_VALUES));
3457:       }
3458:     }
3459:   } else {
3460:     PetscCall(MatSeqAIJGetArrayWrite(x, &aa));
3461:     for (i = 0; i < aij->nz; i++) PetscCall(PetscRandomGetValue(rctx, aa + i));
3462:     PetscCall(MatSeqAIJRestoreArrayWrite(x, &aa));
3463:   }
3464:   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
3465:   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
3466:   PetscFunctionReturn(PETSC_SUCCESS);
3467: }

3469: /* Like MatSetRandom_SeqAIJ, but do not set values on columns in range of [low, high) */
3470: PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat x, PetscInt low, PetscInt high, PetscRandom rctx)
3471: {
3472:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)x->data;
3473:   PetscScalar a;
3474:   PetscInt    m, n, i, j, col, nskip;

3476:   PetscFunctionBegin;
3477:   nskip = high - low;
3478:   PetscCall(MatGetSize(x, &m, &n));
3479:   n -= nskip; /* shrink number of columns where nonzeros can be set */
3480:   for (i = 0; i < m; i++) {
3481:     for (j = 0; j < aij->imax[i]; j++) {
3482:       PetscCall(PetscRandomGetValue(rctx, &a));
3483:       col = (PetscInt)(n * PetscRealPart(a));
3484:       if (col >= low) col += nskip; /* shift col rightward to skip the hole */
3485:       PetscCall(MatSetValues(x, 1, &i, 1, &col, &a, ADD_VALUES));
3486:     }
3487:   }
3488:   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
3489:   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
3490:   PetscFunctionReturn(PETSC_SUCCESS);
3491: }

3493: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
3494:                                        MatGetRow_SeqAIJ,
3495:                                        MatRestoreRow_SeqAIJ,
3496:                                        MatMult_SeqAIJ,
3497:                                        /*  4*/ MatMultAdd_SeqAIJ,
3498:                                        MatMultTranspose_SeqAIJ,
3499:                                        MatMultTransposeAdd_SeqAIJ,
3500:                                        NULL,
3501:                                        NULL,
3502:                                        NULL,
3503:                                        /* 10*/ NULL,
3504:                                        MatLUFactor_SeqAIJ,
3505:                                        NULL,
3506:                                        MatSOR_SeqAIJ,
3507:                                        MatTranspose_SeqAIJ,
3508:                                        /*1 5*/ MatGetInfo_SeqAIJ,
3509:                                        MatEqual_SeqAIJ,
3510:                                        MatGetDiagonal_SeqAIJ,
3511:                                        MatDiagonalScale_SeqAIJ,
3512:                                        MatNorm_SeqAIJ,
3513:                                        /* 20*/ NULL,
3514:                                        MatAssemblyEnd_SeqAIJ,
3515:                                        MatSetOption_SeqAIJ,
3516:                                        MatZeroEntries_SeqAIJ,
3517:                                        /* 24*/ MatZeroRows_SeqAIJ,
3518:                                        NULL,
3519:                                        NULL,
3520:                                        NULL,
3521:                                        NULL,
3522:                                        /* 29*/ MatSetUp_Seq_Hash,
3523:                                        NULL,
3524:                                        NULL,
3525:                                        NULL,
3526:                                        NULL,
3527:                                        /* 34*/ MatDuplicate_SeqAIJ,
3528:                                        NULL,
3529:                                        NULL,
3530:                                        MatILUFactor_SeqAIJ,
3531:                                        NULL,
3532:                                        /* 39*/ MatAXPY_SeqAIJ,
3533:                                        MatCreateSubMatrices_SeqAIJ,
3534:                                        MatIncreaseOverlap_SeqAIJ,
3535:                                        MatGetValues_SeqAIJ,
3536:                                        MatCopy_SeqAIJ,
3537:                                        /* 44*/ MatGetRowMax_SeqAIJ,
3538:                                        MatScale_SeqAIJ,
3539:                                        MatShift_SeqAIJ,
3540:                                        MatDiagonalSet_SeqAIJ,
3541:                                        MatZeroRowsColumns_SeqAIJ,
3542:                                        /* 49*/ MatSetRandom_SeqAIJ,
3543:                                        MatGetRowIJ_SeqAIJ,
3544:                                        MatRestoreRowIJ_SeqAIJ,
3545:                                        MatGetColumnIJ_SeqAIJ,
3546:                                        MatRestoreColumnIJ_SeqAIJ,
3547:                                        /* 54*/ MatFDColoringCreate_SeqXAIJ,
3548:                                        NULL,
3549:                                        NULL,
3550:                                        MatPermute_SeqAIJ,
3551:                                        NULL,
3552:                                        /* 59*/ NULL,
3553:                                        MatDestroy_SeqAIJ,
3554:                                        MatView_SeqAIJ,
3555:                                        NULL,
3556:                                        NULL,
3557:                                        /* 64*/ NULL,
3558:                                        MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ,
3559:                                        NULL,
3560:                                        NULL,
3561:                                        NULL,
3562:                                        /* 69*/ MatGetRowMaxAbs_SeqAIJ,
3563:                                        MatGetRowMinAbs_SeqAIJ,
3564:                                        NULL,
3565:                                        NULL,
3566:                                        NULL,
3567:                                        /* 74*/ NULL,
3568:                                        MatFDColoringApply_AIJ,
3569:                                        NULL,
3570:                                        NULL,
3571:                                        NULL,
3572:                                        /* 79*/ MatFindZeroDiagonals_SeqAIJ,
3573:                                        NULL,
3574:                                        NULL,
3575:                                        NULL,
3576:                                        MatLoad_SeqAIJ,
3577:                                        /* 84*/ NULL,
3578:                                        NULL,
3579:                                        NULL,
3580:                                        NULL,
3581:                                        NULL,
3582:                                        /* 89*/ NULL,
3583:                                        NULL,
3584:                                        MatMatMultNumeric_SeqAIJ_SeqAIJ,
3585:                                        NULL,
3586:                                        NULL,
3587:                                        /* 94*/ MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy,
3588:                                        NULL,
3589:                                        NULL,
3590:                                        MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ,
3591:                                        NULL,
3592:                                        /* 99*/ MatProductSetFromOptions_SeqAIJ,
3593:                                        NULL,
3594:                                        NULL,
3595:                                        MatConjugate_SeqAIJ,
3596:                                        NULL,
3597:                                        /*104*/ MatSetValuesRow_SeqAIJ,
3598:                                        MatRealPart_SeqAIJ,
3599:                                        MatImaginaryPart_SeqAIJ,
3600:                                        NULL,
3601:                                        NULL,
3602:                                        /*109*/ MatMatSolve_SeqAIJ,
3603:                                        NULL,
3604:                                        MatGetRowMin_SeqAIJ,
3605:                                        NULL,
3606:                                        MatMissingDiagonal_SeqAIJ,
3607:                                        /*114*/ NULL,
3608:                                        NULL,
3609:                                        NULL,
3610:                                        NULL,
3611:                                        NULL,
3612:                                        /*119*/ NULL,
3613:                                        NULL,
3614:                                        NULL,
3615:                                        NULL,
3616:                                        MatGetMultiProcBlock_SeqAIJ,
3617:                                        /*124*/ MatFindNonzeroRows_SeqAIJ,
3618:                                        MatGetColumnReductions_SeqAIJ,
3619:                                        MatInvertBlockDiagonal_SeqAIJ,
3620:                                        MatInvertVariableBlockDiagonal_SeqAIJ,
3621:                                        NULL,
3622:                                        /*129*/ NULL,
3623:                                        NULL,
3624:                                        NULL,
3625:                                        MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ,
3626:                                        MatTransposeColoringCreate_SeqAIJ,
3627:                                        /*134*/ MatTransColoringApplySpToDen_SeqAIJ,
3628:                                        MatTransColoringApplyDenToSp_SeqAIJ,
3629:                                        NULL,
3630:                                        NULL,
3631:                                        MatRARtNumeric_SeqAIJ_SeqAIJ,
3632:                                        /*139*/ NULL,
3633:                                        NULL,
3634:                                        NULL,
3635:                                        MatFDColoringSetUp_SeqXAIJ,
3636:                                        MatFindOffBlockDiagonalEntries_SeqAIJ,
3637:                                        MatCreateMPIMatConcatenateSeqMat_SeqAIJ,
3638:                                        /*145*/ MatDestroySubMatrices_SeqAIJ,
3639:                                        NULL,
3640:                                        NULL,
3641:                                        MatCreateGraph_Simple_AIJ,
3642:                                        NULL,
3643:                                        /*150*/ MatTransposeSymbolic_SeqAIJ,
3644:                                        MatEliminateZeros_SeqAIJ,
3645:                                        MatGetRowSumAbs_SeqAIJ,
3646:                                        NULL,
3647:                                        NULL,
3648:                                        /*155*/ NULL,
3649:                                        MatCopyHashToXAIJ_Seq_Hash};

3651: static PetscErrorCode MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat, PetscInt *indices)
3652: {
3653:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
3654:   PetscInt    i, nz, n;

3656:   PetscFunctionBegin;
3657:   nz = aij->maxnz;
3658:   n  = mat->rmap->n;
3659:   for (i = 0; i < nz; i++) aij->j[i] = indices[i];
3660:   aij->nz = nz;
3661:   for (i = 0; i < n; i++) aij->ilen[i] = aij->imax[i];
3662:   PetscFunctionReturn(PETSC_SUCCESS);
3663: }

3665: /*
3666:  * Given a sparse matrix with global column indices, compact it by using a local column space.
3667:  * The result matrix helps saving memory in other algorithms, such as MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable()
3668:  */
3669: PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat mat, ISLocalToGlobalMapping *mapping)
3670: {
3671:   Mat_SeqAIJ   *aij = (Mat_SeqAIJ *)mat->data;
3672:   PetscHMapI    gid1_lid1;
3673:   PetscHashIter tpos;
3674:   PetscInt      gid, lid, i, ec, nz = aij->nz;
3675:   PetscInt     *garray, *jj = aij->j;

3677:   PetscFunctionBegin;
3679:   PetscAssertPointer(mapping, 2);
3680:   /* use a table */
3681:   PetscCall(PetscHMapICreateWithSize(mat->rmap->n, &gid1_lid1));
3682:   ec = 0;
3683:   for (i = 0; i < nz; i++) {
3684:     PetscInt data, gid1 = jj[i] + 1;
3685:     PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &data));
3686:     if (!data) {
3687:       /* one based table */
3688:       PetscCall(PetscHMapISet(gid1_lid1, gid1, ++ec));
3689:     }
3690:   }
3691:   /* form array of columns we need */
3692:   PetscCall(PetscMalloc1(ec, &garray));
3693:   PetscHashIterBegin(gid1_lid1, tpos);
3694:   while (!PetscHashIterAtEnd(gid1_lid1, tpos)) {
3695:     PetscHashIterGetKey(gid1_lid1, tpos, gid);
3696:     PetscHashIterGetVal(gid1_lid1, tpos, lid);
3697:     PetscHashIterNext(gid1_lid1, tpos);
3698:     gid--;
3699:     lid--;
3700:     garray[lid] = gid;
3701:   }
3702:   PetscCall(PetscSortInt(ec, garray)); /* sort, and rebuild */
3703:   PetscCall(PetscHMapIClear(gid1_lid1));
3704:   for (i = 0; i < ec; i++) PetscCall(PetscHMapISet(gid1_lid1, garray[i] + 1, i + 1));
3705:   /* compact out the extra columns in B */
3706:   for (i = 0; i < nz; i++) {
3707:     PetscInt gid1 = jj[i] + 1;
3708:     PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &lid));
3709:     lid--;
3710:     jj[i] = lid;
3711:   }
3712:   PetscCall(PetscLayoutDestroy(&mat->cmap));
3713:   PetscCall(PetscHMapIDestroy(&gid1_lid1));
3714:   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)mat), ec, ec, 1, &mat->cmap));
3715:   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, mat->cmap->bs, mat->cmap->n, garray, PETSC_OWN_POINTER, mapping));
3716:   PetscCall(ISLocalToGlobalMappingSetType(*mapping, ISLOCALTOGLOBALMAPPINGHASH));
3717:   PetscFunctionReturn(PETSC_SUCCESS);
3718: }

3720: /*@
3721:   MatSeqAIJSetColumnIndices - Set the column indices for all the rows
3722:   in the matrix.

3724:   Input Parameters:
3725: + mat     - the `MATSEQAIJ` matrix
3726: - indices - the column indices

3728:   Level: advanced

3730:   Notes:
3731:   This can be called if you have precomputed the nonzero structure of the
3732:   matrix and want to provide it to the matrix object to improve the performance
3733:   of the `MatSetValues()` operation.

3735:   You MUST have set the correct numbers of nonzeros per row in the call to
3736:   `MatCreateSeqAIJ()`, and the columns indices MUST be sorted.

3738:   MUST be called before any calls to `MatSetValues()`

3740:   The indices should start with zero, not one.

3742: .seealso: [](ch_matrices), `Mat`, `MATSEQAIJ`
3743: @*/
3744: PetscErrorCode MatSeqAIJSetColumnIndices(Mat mat, PetscInt *indices)
3745: {
3746:   PetscFunctionBegin;
3748:   PetscAssertPointer(indices, 2);
3749:   PetscUseMethod(mat, "MatSeqAIJSetColumnIndices_C", (Mat, PetscInt *), (mat, indices));
3750:   PetscFunctionReturn(PETSC_SUCCESS);
3751: }

3753: static PetscErrorCode MatStoreValues_SeqAIJ(Mat mat)
3754: {
3755:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
3756:   size_t      nz  = aij->i[mat->rmap->n];

3758:   PetscFunctionBegin;
3759:   PetscCheck(aij->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");

3761:   /* allocate space for values if not already there */
3762:   if (!aij->saved_values) { PetscCall(PetscMalloc1(nz + 1, &aij->saved_values)); }

3764:   /* copy values over */
3765:   PetscCall(PetscArraycpy(aij->saved_values, aij->a, nz));
3766:   PetscFunctionReturn(PETSC_SUCCESS);
3767: }

3769: /*@
3770:   MatStoreValues - Stashes a copy of the matrix values; this allows reusing of the linear part of a Jacobian, while recomputing only the
3771:   nonlinear portion.

3773:   Logically Collect

3775:   Input Parameter:
3776: . mat - the matrix (currently only `MATAIJ` matrices support this option)

3778:   Level: advanced

3780:   Example Usage:
3781: .vb
3782:     Using SNES
3783:     Create Jacobian matrix
3784:     Set linear terms into matrix
3785:     Apply boundary conditions to matrix, at this time matrix must have
3786:       final nonzero structure (i.e. setting the nonlinear terms and applying
3787:       boundary conditions again will not change the nonzero structure
3788:     MatSetOption(mat, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE);
3789:     MatStoreValues(mat);
3790:     Call SNESSetJacobian() with matrix
3791:     In your Jacobian routine
3792:       MatRetrieveValues(mat);
3793:       Set nonlinear terms in matrix

3795:     Without `SNESSolve()`, i.e. when you handle nonlinear solve yourself:
3796:     // build linear portion of Jacobian
3797:     MatSetOption(mat, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE);
3798:     MatStoreValues(mat);
3799:     loop over nonlinear iterations
3800:        MatRetrieveValues(mat);
3801:        // call MatSetValues(mat,...) to set nonliner portion of Jacobian
3802:        // call MatAssemblyBegin/End() on matrix
3803:        Solve linear system with Jacobian
3804:     endloop
3805: .ve

3807:   Notes:
3808:   Matrix must already be assembled before calling this routine
3809:   Must set the matrix option `MatSetOption`(mat,`MAT_NEW_NONZERO_LOCATIONS`,`PETSC_FALSE`); before
3810:   calling this routine.

3812:   When this is called multiple times it overwrites the previous set of stored values
3813:   and does not allocated additional space.

3815: .seealso: [](ch_matrices), `Mat`, `MatRetrieveValues()`
3816: @*/
3817: PetscErrorCode MatStoreValues(Mat mat)
3818: {
3819:   PetscFunctionBegin;
3821:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3822:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3823:   PetscUseMethod(mat, "MatStoreValues_C", (Mat), (mat));
3824:   PetscFunctionReturn(PETSC_SUCCESS);
3825: }

3827: static PetscErrorCode MatRetrieveValues_SeqAIJ(Mat mat)
3828: {
3829:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
3830:   PetscInt    nz  = aij->i[mat->rmap->n];

3832:   PetscFunctionBegin;
3833:   PetscCheck(aij->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
3834:   PetscCheck(aij->saved_values, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatStoreValues(A);first");
3835:   /* copy values over */
3836:   PetscCall(PetscArraycpy(aij->a, aij->saved_values, nz));
3837:   PetscFunctionReturn(PETSC_SUCCESS);
3838: }

3840: /*@
3841:   MatRetrieveValues - Retrieves the copy of the matrix values that was stored with `MatStoreValues()`

3843:   Logically Collect

3845:   Input Parameter:
3846: . mat - the matrix (currently only `MATAIJ` matrices support this option)

3848:   Level: advanced

3850: .seealso: [](ch_matrices), `Mat`, `MatStoreValues()`
3851: @*/
3852: PetscErrorCode MatRetrieveValues(Mat mat)
3853: {
3854:   PetscFunctionBegin;
3856:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3857:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3858:   PetscUseMethod(mat, "MatRetrieveValues_C", (Mat), (mat));
3859:   PetscFunctionReturn(PETSC_SUCCESS);
3860: }

3862: /*@
3863:   MatCreateSeqAIJ - Creates a sparse matrix in `MATSEQAIJ` (compressed row) format
3864:   (the default parallel PETSc format).  For good matrix assembly performance
3865:   the user should preallocate the matrix storage by setting the parameter `nz`
3866:   (or the array `nnz`).

3868:   Collective

3870:   Input Parameters:
3871: + comm - MPI communicator, set to `PETSC_COMM_SELF`
3872: . m    - number of rows
3873: . n    - number of columns
3874: . nz   - number of nonzeros per row (same for all rows)
3875: - nnz  - array containing the number of nonzeros in the various rows
3876:          (possibly different for each row) or NULL

3878:   Output Parameter:
3879: . A - the matrix

3881:   Options Database Keys:
3882: + -mat_no_inode            - Do not use inodes
3883: - -mat_inode_limit <limit> - Sets inode limit (max limit=5)

3885:   Level: intermediate

3887:   Notes:
3888:   It is recommend to use `MatCreateFromOptions()` instead of this routine

3890:   If `nnz` is given then `nz` is ignored

3892:   The `MATSEQAIJ` format, also called
3893:   compressed row storage, is fully compatible with standard Fortran
3894:   storage.  That is, the stored row and column indices can begin at
3895:   either one (as in Fortran) or zero.

3897:   Specify the preallocated storage with either `nz` or `nnz` (not both).
3898:   Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
3899:   allocation.

3901:   By default, this format uses inodes (identical nodes) when possible, to
3902:   improve numerical efficiency of matrix-vector products and solves. We
3903:   search for consecutive rows with the same nonzero structure, thereby
3904:   reusing matrix information to achieve increased efficiency.

3906: .seealso: [](ch_matrices), `Mat`, [Sparse Matrix Creation](sec_matsparse), `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()`
3907: @*/
3908: PetscErrorCode MatCreateSeqAIJ(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A)
3909: {
3910:   PetscFunctionBegin;
3911:   PetscCall(MatCreate(comm, A));
3912:   PetscCall(MatSetSizes(*A, m, n, m, n));
3913:   PetscCall(MatSetType(*A, MATSEQAIJ));
3914:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*A, nz, nnz));
3915:   PetscFunctionReturn(PETSC_SUCCESS);
3916: }

3918: /*@
3919:   MatSeqAIJSetPreallocation - For good matrix assembly performance
3920:   the user should preallocate the matrix storage by setting the parameter nz
3921:   (or the array nnz).  By setting these parameters accurately, performance
3922:   during matrix assembly can be increased by more than a factor of 50.

3924:   Collective

3926:   Input Parameters:
3927: + B   - The matrix
3928: . nz  - number of nonzeros per row (same for all rows)
3929: - nnz - array containing the number of nonzeros in the various rows
3930:          (possibly different for each row) or NULL

3932:   Options Database Keys:
3933: + -mat_no_inode            - Do not use inodes
3934: - -mat_inode_limit <limit> - Sets inode limit (max limit=5)

3936:   Level: intermediate

3938:   Notes:
3939:   If `nnz` is given then `nz` is ignored

3941:   The `MATSEQAIJ` format also called
3942:   compressed row storage, is fully compatible with standard Fortran
3943:   storage.  That is, the stored row and column indices can begin at
3944:   either one (as in Fortran) or zero.  See the users' manual for details.

3946:   Specify the preallocated storage with either `nz` or `nnz` (not both).
3947:   Set nz = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
3948:   allocation.

3950:   You can call `MatGetInfo()` to get information on how effective the preallocation was;
3951:   for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
3952:   You can also run with the option -info and look for messages with the string
3953:   malloc in them to see if additional memory allocation was needed.

3955:   Developer Notes:
3956:   Use nz of `MAT_SKIP_ALLOCATION` to not allocate any space for the matrix
3957:   entries or columns indices

3959:   By default, this format uses inodes (identical nodes) when possible, to
3960:   improve numerical efficiency of matrix-vector products and solves. We
3961:   search for consecutive rows with the same nonzero structure, thereby
3962:   reusing matrix information to achieve increased efficiency.

3964: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()`, `MatGetInfo()`,
3965:           `MatSeqAIJSetTotalPreallocation()`
3966: @*/
3967: PetscErrorCode MatSeqAIJSetPreallocation(Mat B, PetscInt nz, const PetscInt nnz[])
3968: {
3969:   PetscFunctionBegin;
3972:   PetscTryMethod(B, "MatSeqAIJSetPreallocation_C", (Mat, PetscInt, const PetscInt[]), (B, nz, nnz));
3973:   PetscFunctionReturn(PETSC_SUCCESS);
3974: }

3976: PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat B, PetscInt nz, const PetscInt *nnz)
3977: {
3978:   Mat_SeqAIJ *b              = (Mat_SeqAIJ *)B->data;
3979:   PetscBool   skipallocation = PETSC_FALSE, realalloc = PETSC_FALSE;
3980:   PetscInt    i;

3982:   PetscFunctionBegin;
3983:   if (B->hash_active) {
3984:     B->ops[0] = b->cops;
3985:     PetscCall(PetscHMapIJVDestroy(&b->ht));
3986:     PetscCall(PetscFree(b->dnz));
3987:     B->hash_active = PETSC_FALSE;
3988:   }
3989:   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
3990:   if (nz == MAT_SKIP_ALLOCATION) {
3991:     skipallocation = PETSC_TRUE;
3992:     nz             = 0;
3993:   }
3994:   PetscCall(PetscLayoutSetUp(B->rmap));
3995:   PetscCall(PetscLayoutSetUp(B->cmap));

3997:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
3998:   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nz cannot be less than 0: value %" PetscInt_FMT, nz);
3999:   if (nnz) {
4000:     for (i = 0; i < B->rmap->n; i++) {
4001:       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]);
4002:       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);
4003:     }
4004:   }

4006:   B->preallocated = PETSC_TRUE;
4007:   if (!skipallocation) {
4008:     if (!b->imax) { PetscCall(PetscMalloc1(B->rmap->n, &b->imax)); }
4009:     if (!b->ilen) {
4010:       /* b->ilen will count nonzeros in each row so far. */
4011:       PetscCall(PetscCalloc1(B->rmap->n, &b->ilen));
4012:     } else {
4013:       PetscCall(PetscMemzero(b->ilen, B->rmap->n * sizeof(PetscInt)));
4014:     }
4015:     if (!b->ipre) PetscCall(PetscMalloc1(B->rmap->n, &b->ipre));
4016:     if (!nnz) {
4017:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
4018:       else if (nz < 0) nz = 1;
4019:       nz = PetscMin(nz, B->cmap->n);
4020:       for (i = 0; i < B->rmap->n; i++) b->imax[i] = nz;
4021:       PetscCall(PetscIntMultError(nz, B->rmap->n, &nz));
4022:     } else {
4023:       PetscInt64 nz64 = 0;
4024:       for (i = 0; i < B->rmap->n; i++) {
4025:         b->imax[i] = nnz[i];
4026:         nz64 += nnz[i];
4027:       }
4028:       PetscCall(PetscIntCast(nz64, &nz));
4029:     }

4031:     /* allocate the matrix space */
4032:     PetscCall(MatSeqXAIJFreeAIJ(B, &b->a, &b->j, &b->i));
4033:     PetscCall(PetscShmgetAllocateArray(nz, sizeof(PetscInt), (void **)&b->j));
4034:     PetscCall(PetscShmgetAllocateArray(B->rmap->n + 1, sizeof(PetscInt), (void **)&b->i));
4035:     b->free_ij = PETSC_TRUE;
4036:     if (B->structure_only) {
4037:       b->free_a = PETSC_FALSE;
4038:     } else {
4039:       PetscCall(PetscShmgetAllocateArray(nz, sizeof(PetscScalar), (void **)&b->a));
4040:       b->free_a = PETSC_TRUE;
4041:     }
4042:     b->i[0] = 0;
4043:     for (i = 1; i < B->rmap->n + 1; i++) b->i[i] = b->i[i - 1] + b->imax[i - 1];
4044:   } else {
4045:     b->free_a  = PETSC_FALSE;
4046:     b->free_ij = PETSC_FALSE;
4047:   }

4049:   if (b->ipre && nnz != b->ipre && b->imax) {
4050:     /* reserve user-requested sparsity */
4051:     PetscCall(PetscArraycpy(b->ipre, b->imax, B->rmap->n));
4052:   }

4054:   b->nz               = 0;
4055:   b->maxnz            = nz;
4056:   B->info.nz_unneeded = (double)b->maxnz;
4057:   if (realalloc) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
4058:   B->was_assembled = PETSC_FALSE;
4059:   B->assembled     = PETSC_FALSE;
4060:   /* We simply deem preallocation has changed nonzero state. Updating the state
4061:      will give clients (like AIJKokkos) a chance to know something has happened.
4062:   */
4063:   B->nonzerostate++;
4064:   PetscFunctionReturn(PETSC_SUCCESS);
4065: }

4067: static PetscErrorCode MatResetPreallocation_SeqAIJ(Mat A)
4068: {
4069:   Mat_SeqAIJ *a;
4070:   PetscInt    i;
4071:   PetscBool   skipreset;

4073:   PetscFunctionBegin;

4076:   /* Check local size. If zero, then return */
4077:   if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);

4079:   a = (Mat_SeqAIJ *)A->data;
4080:   /* if no saved info, we error out */
4081:   PetscCheck(a->ipre, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "No saved preallocation info ");

4083:   PetscCheck(a->i && a->imax && a->ilen, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Memory info is incomplete, and can not reset preallocation ");

4085:   PetscCall(PetscArraycmp(a->ipre, a->ilen, A->rmap->n, &skipreset));
4086:   if (!skipreset) {
4087:     PetscCall(PetscArraycpy(a->imax, a->ipre, A->rmap->n));
4088:     PetscCall(PetscArrayzero(a->ilen, A->rmap->n));
4089:     a->i[0] = 0;
4090:     for (i = 1; i < A->rmap->n + 1; i++) a->i[i] = a->i[i - 1] + a->imax[i - 1];
4091:     A->preallocated     = PETSC_TRUE;
4092:     a->nz               = 0;
4093:     a->maxnz            = a->i[A->rmap->n];
4094:     A->info.nz_unneeded = (double)a->maxnz;
4095:     A->was_assembled    = PETSC_FALSE;
4096:     A->assembled        = PETSC_FALSE;
4097:   }
4098:   PetscFunctionReturn(PETSC_SUCCESS);
4099: }

4101: /*@
4102:   MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in `MATSEQAIJ` format.

4104:   Input Parameters:
4105: + B - the matrix
4106: . i - the indices into `j` for the start of each row (indices start with zero)
4107: . j - the column indices for each row (indices start with zero) these must be sorted for each row
4108: - v - optional values in the matrix, use `NULL` if not provided

4110:   Level: developer

4112:   Notes:
4113:   The `i`,`j`,`v` values are COPIED with this routine; to avoid the copy use `MatCreateSeqAIJWithArrays()`

4115:   This routine may be called multiple times with different nonzero patterns (or the same nonzero pattern). The nonzero
4116:   structure will be the union of all the previous nonzero structures.

4118:   Developer Notes:
4119:   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
4120:   then just copies the `v` values directly with `PetscMemcpy()`.

4122:   This routine could also take a `PetscCopyMode` argument to allow sharing the values instead of always copying them.

4124: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatSeqAIJSetPreallocation()`, `MATSEQAIJ`, `MatResetPreallocation()`
4125: @*/
4126: PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B, const PetscInt i[], const PetscInt j[], const PetscScalar v[])
4127: {
4128:   PetscFunctionBegin;
4131:   PetscTryMethod(B, "MatSeqAIJSetPreallocationCSR_C", (Mat, const PetscInt[], const PetscInt[], const PetscScalar[]), (B, i, j, v));
4132:   PetscFunctionReturn(PETSC_SUCCESS);
4133: }

4135: static PetscErrorCode MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B, const PetscInt Ii[], const PetscInt J[], const PetscScalar v[])
4136: {
4137:   PetscInt  i;
4138:   PetscInt  m, n;
4139:   PetscInt  nz;
4140:   PetscInt *nnz;

4142:   PetscFunctionBegin;
4143:   PetscCheck(Ii[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Ii[0] must be 0 it is %" PetscInt_FMT, Ii[0]);

4145:   PetscCall(PetscLayoutSetUp(B->rmap));
4146:   PetscCall(PetscLayoutSetUp(B->cmap));

4148:   PetscCall(MatGetSize(B, &m, &n));
4149:   PetscCall(PetscMalloc1(m + 1, &nnz));
4150:   for (i = 0; i < m; i++) {
4151:     nz = Ii[i + 1] - Ii[i];
4152:     PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local row %" PetscInt_FMT " has a negative number of columns %" PetscInt_FMT, i, nz);
4153:     nnz[i] = nz;
4154:   }
4155:   PetscCall(MatSeqAIJSetPreallocation(B, 0, nnz));
4156:   PetscCall(PetscFree(nnz));

4158:   for (i = 0; i < m; i++) PetscCall(MatSetValues_SeqAIJ(B, 1, &i, Ii[i + 1] - Ii[i], J + Ii[i], PetscSafePointerPlusOffset(v, Ii[i]), INSERT_VALUES));

4160:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
4161:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));

4163:   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
4164:   PetscFunctionReturn(PETSC_SUCCESS);
4165: }

4167: /*@
4168:   MatSeqAIJKron - Computes `C`, the Kronecker product of `A` and `B`.

4170:   Input Parameters:
4171: + A     - left-hand side matrix
4172: . B     - right-hand side matrix
4173: - reuse - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`

4175:   Output Parameter:
4176: . C - Kronecker product of `A` and `B`

4178:   Level: intermediate

4180:   Note:
4181:   `MAT_REUSE_MATRIX` can only be used when the nonzero structure of the product matrix has not changed from that last call to `MatSeqAIJKron()`.

4183: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MATSEQAIJ`, `MATKAIJ`, `MatReuse`
4184: @*/
4185: PetscErrorCode MatSeqAIJKron(Mat A, Mat B, MatReuse reuse, Mat *C)
4186: {
4187:   PetscFunctionBegin;
4192:   PetscAssertPointer(C, 4);
4193:   if (reuse == MAT_REUSE_MATRIX) {
4196:   }
4197:   PetscTryMethod(A, "MatSeqAIJKron_C", (Mat, Mat, MatReuse, Mat *), (A, B, reuse, C));
4198:   PetscFunctionReturn(PETSC_SUCCESS);
4199: }

4201: static PetscErrorCode MatSeqAIJKron_SeqAIJ(Mat A, Mat B, MatReuse reuse, Mat *C)
4202: {
4203:   Mat                newmat;
4204:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
4205:   Mat_SeqAIJ        *b = (Mat_SeqAIJ *)B->data;
4206:   PetscScalar       *v;
4207:   const PetscScalar *aa, *ba;
4208:   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;
4209:   PetscBool          flg;

4211:   PetscFunctionBegin;
4212:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4213:   PetscCheck(A->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
4214:   PetscCheck(!B->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4215:   PetscCheck(B->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
4216:   PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &flg));
4217:   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatType %s", ((PetscObject)B)->type_name);
4218:   PetscCheck(reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatReuse %d", (int)reuse);
4219:   if (reuse == MAT_INITIAL_MATRIX) {
4220:     PetscCall(PetscMalloc2(am * bm + 1, &i, a->i[am] * b->i[bm], &j));
4221:     PetscCall(MatCreate(PETSC_COMM_SELF, &newmat));
4222:     PetscCall(MatSetSizes(newmat, am * bm, an * bn, am * bm, an * bn));
4223:     PetscCall(MatSetType(newmat, MATAIJ));
4224:     i[0] = 0;
4225:     for (m = 0; m < am; ++m) {
4226:       for (p = 0; p < bm; ++p) {
4227:         i[m * bm + p + 1] = i[m * bm + p] + (a->i[m + 1] - a->i[m]) * (b->i[p + 1] - b->i[p]);
4228:         for (n = a->i[m]; n < a->i[m + 1]; ++n) {
4229:           for (q = b->i[p]; q < b->i[p + 1]; ++q) j[nnz++] = a->j[n] * bn + b->j[q];
4230:         }
4231:       }
4232:     }
4233:     PetscCall(MatSeqAIJSetPreallocationCSR(newmat, i, j, NULL));
4234:     *C = newmat;
4235:     PetscCall(PetscFree2(i, j));
4236:     nnz = 0;
4237:   }
4238:   PetscCall(MatSeqAIJGetArray(*C, &v));
4239:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
4240:   PetscCall(MatSeqAIJGetArrayRead(B, &ba));
4241:   for (m = 0; m < am; ++m) {
4242:     for (p = 0; p < bm; ++p) {
4243:       for (n = a->i[m]; n < a->i[m + 1]; ++n) {
4244:         for (q = b->i[p]; q < b->i[p + 1]; ++q) v[nnz++] = aa[n] * ba[q];
4245:       }
4246:     }
4247:   }
4248:   PetscCall(MatSeqAIJRestoreArray(*C, &v));
4249:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
4250:   PetscCall(MatSeqAIJRestoreArrayRead(B, &ba));
4251:   PetscFunctionReturn(PETSC_SUCCESS);
4252: }

4254: #include <../src/mat/impls/dense/seq/dense.h>
4255: #include <petsc/private/kernels/petscaxpy.h>

4257: /*
4258:     Computes (B'*A')' since computing B*A directly is untenable

4260:                n                       p                          p
4261:         [             ]       [             ]         [                 ]
4262:       m [      A      ]  *  n [       B     ]   =   m [         C       ]
4263:         [             ]       [             ]         [                 ]

4265: */
4266: PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat A, Mat B, Mat C)
4267: {
4268:   Mat_SeqDense      *sub_a = (Mat_SeqDense *)A->data;
4269:   Mat_SeqAIJ        *sub_b = (Mat_SeqAIJ *)B->data;
4270:   Mat_SeqDense      *sub_c = (Mat_SeqDense *)C->data;
4271:   PetscInt           i, j, n, m, q, p;
4272:   const PetscInt    *ii, *idx;
4273:   const PetscScalar *b, *a, *a_q;
4274:   PetscScalar       *c, *c_q;
4275:   PetscInt           clda = sub_c->lda;
4276:   PetscInt           alda = sub_a->lda;

4278:   PetscFunctionBegin;
4279:   m = A->rmap->n;
4280:   n = A->cmap->n;
4281:   p = B->cmap->n;
4282:   a = sub_a->v;
4283:   b = sub_b->a;
4284:   c = sub_c->v;
4285:   if (clda == m) {
4286:     PetscCall(PetscArrayzero(c, m * p));
4287:   } else {
4288:     for (j = 0; j < p; j++)
4289:       for (i = 0; i < m; i++) c[j * clda + i] = 0.0;
4290:   }
4291:   ii  = sub_b->i;
4292:   idx = sub_b->j;
4293:   for (i = 0; i < n; i++) {
4294:     q = ii[i + 1] - ii[i];
4295:     while (q-- > 0) {
4296:       c_q = c + clda * (*idx);
4297:       a_q = a + alda * i;
4298:       PetscKernelAXPY(c_q, *b, a_q, m);
4299:       idx++;
4300:       b++;
4301:     }
4302:   }
4303:   PetscFunctionReturn(PETSC_SUCCESS);
4304: }

4306: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat A, Mat B, PetscReal fill, Mat C)
4307: {
4308:   PetscInt  m = A->rmap->n, n = B->cmap->n;
4309:   PetscBool cisdense;

4311:   PetscFunctionBegin;
4312:   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);
4313:   PetscCall(MatSetSizes(C, m, n, m, n));
4314:   PetscCall(MatSetBlockSizesFromMats(C, A, B));
4315:   PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATSEQDENSECUDA, MATSEQDENSEHIP, ""));
4316:   if (!cisdense) PetscCall(MatSetType(C, MATDENSE));
4317:   PetscCall(MatSetUp(C));

4319:   C->ops->matmultnumeric = MatMatMultNumeric_SeqDense_SeqAIJ;
4320:   PetscFunctionReturn(PETSC_SUCCESS);
4321: }

4323: /*MC
4324:    MATSEQAIJ - MATSEQAIJ = "seqaij" - A matrix type to be used for sequential sparse matrices,
4325:    based on compressed sparse row format.

4327:    Options Database Key:
4328: . -mat_type seqaij - sets the matrix type to "seqaij" during a call to MatSetFromOptions()

4330:    Level: beginner

4332:    Notes:
4333:     `MatSetValues()` may be called for this matrix type with a `NULL` argument for the numerical values,
4334:     in this case the values associated with the rows and columns one passes in are set to zero
4335:     in the matrix

4337:     `MatSetOptions`(,`MAT_STRUCTURE_ONLY`,`PETSC_TRUE`) may be called for this matrix type. In this no
4338:     space is allocated for the nonzero entries and any entries passed with `MatSetValues()` are ignored

4340:   Developer Note:
4341:     It would be nice if all matrix formats supported passing `NULL` in for the numerical values

4343: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatSetFromOptions()`, `MatSetType()`, `MatCreate()`, `MatType`, `MATSELL`, `MATSEQSELL`, `MATMPISELL`
4344: M*/

4346: /*MC
4347:    MATAIJ - MATAIJ = "aij" - A matrix type to be used for sparse matrices.

4349:    This matrix type is identical to `MATSEQAIJ` when constructed with a single process communicator,
4350:    and `MATMPIAIJ` otherwise.  As a result, for single process communicators,
4351:    `MatSeqAIJSetPreallocation()` is supported, and similarly `MatMPIAIJSetPreallocation()` is supported
4352:    for communicators controlling multiple processes.  It is recommended that you call both of
4353:    the above preallocation routines for simplicity.

4355:    Options Database Key:
4356: . -mat_type aij - sets the matrix type to "aij" during a call to `MatSetFromOptions()`

4358:   Level: beginner

4360:    Note:
4361:    Subclasses include `MATAIJCUSPARSE`, `MATAIJPERM`, `MATAIJSELL`, `MATAIJMKL`, `MATAIJCRL`, and also automatically switches over to use inodes when
4362:    enough exist.

4364: .seealso: [](ch_matrices), `Mat`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MATSEQAIJ`, `MATMPIAIJ`, `MATSELL`, `MATSEQSELL`, `MATMPISELL`
4365: M*/

4367: /*MC
4368:    MATAIJCRL - MATAIJCRL = "aijcrl" - A matrix type to be used for sparse matrices.

4370:    Options Database Key:
4371: . -mat_type aijcrl - sets the matrix type to "aijcrl" during a call to `MatSetFromOptions()`

4373:   Level: beginner

4375:    Note:
4376:    This matrix type is identical to `MATSEQAIJCRL` when constructed with a single process communicator,
4377:    and `MATMPIAIJCRL` otherwise.  As a result, for single process communicators,
4378:    `MatSeqAIJSetPreallocation()` is supported, and similarly `MatMPIAIJSetPreallocation()` is supported
4379:    for communicators controlling multiple processes.  It is recommended that you call both of
4380:    the above preallocation routines for simplicity.

4382: .seealso: [](ch_matrices), `Mat`, `MatCreateMPIAIJCRL`, `MATSEQAIJCRL`, `MATMPIAIJCRL`, `MATSEQAIJCRL`, `MATMPIAIJCRL`
4383: M*/

4385: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *);
4386: #if defined(PETSC_HAVE_ELEMENTAL)
4387: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
4388: #endif
4389: #if defined(PETSC_HAVE_SCALAPACK)
4390: PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
4391: #endif
4392: #if defined(PETSC_HAVE_HYPRE)
4393: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType, MatReuse, Mat *);
4394: #endif

4396: PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat, MatType, MatReuse, Mat *);
4397: PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat, MatType, MatReuse, Mat *);
4398: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_IS_XAIJ(Mat);

4400: /*@C
4401:   MatSeqAIJGetArray - gives read/write access to the array where the data for a `MATSEQAIJ` matrix is stored

4403:   Not Collective

4405:   Input Parameter:
4406: . A - a `MATSEQAIJ` matrix

4408:   Output Parameter:
4409: . array - pointer to the data

4411:   Level: intermediate

4413:   Fortran Notes:
4414:   `MatSeqAIJGetArray()` Fortran binding is deprecated (since PETSc 3.19), use `MatSeqAIJGetArrayF90()`

4416: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRestoreArray()`, `MatSeqAIJGetArrayF90()`
4417: @*/
4418: PetscErrorCode MatSeqAIJGetArray(Mat A, PetscScalar *array[])
4419: {
4420:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4422:   PetscFunctionBegin;
4423:   if (aij->ops->getarray) {
4424:     PetscCall((*aij->ops->getarray)(A, array));
4425:   } else {
4426:     *array = aij->a;
4427:   }
4428:   PetscFunctionReturn(PETSC_SUCCESS);
4429: }

4431: /*@C
4432:   MatSeqAIJRestoreArray - returns access to the array where the data for a `MATSEQAIJ` matrix is stored obtained by `MatSeqAIJGetArray()`

4434:   Not Collective

4436:   Input Parameters:
4437: + A     - a `MATSEQAIJ` matrix
4438: - array - pointer to the data

4440:   Level: intermediate

4442:   Fortran Notes:
4443:   `MatSeqAIJRestoreArray()` Fortran binding is deprecated (since PETSc 3.19), use `MatSeqAIJRestoreArrayF90()`

4445: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayF90()`
4446: @*/
4447: PetscErrorCode MatSeqAIJRestoreArray(Mat A, PetscScalar *array[])
4448: {
4449:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4451:   PetscFunctionBegin;
4452:   if (aij->ops->restorearray) {
4453:     PetscCall((*aij->ops->restorearray)(A, array));
4454:   } else {
4455:     *array = NULL;
4456:   }
4457:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
4458:   PetscCall(PetscObjectStateIncrease((PetscObject)A));
4459:   PetscFunctionReturn(PETSC_SUCCESS);
4460: }

4462: /*@C
4463:   MatSeqAIJGetArrayRead - gives read-only access to the array where the data for a `MATSEQAIJ` matrix is stored

4465:   Not Collective; No Fortran Support

4467:   Input Parameter:
4468: . A - a `MATSEQAIJ` matrix

4470:   Output Parameter:
4471: . array - pointer to the data

4473:   Level: intermediate

4475: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayRead()`
4476: @*/
4477: PetscErrorCode MatSeqAIJGetArrayRead(Mat A, const PetscScalar *array[])
4478: {
4479:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4481:   PetscFunctionBegin;
4482:   if (aij->ops->getarrayread) {
4483:     PetscCall((*aij->ops->getarrayread)(A, array));
4484:   } else {
4485:     *array = aij->a;
4486:   }
4487:   PetscFunctionReturn(PETSC_SUCCESS);
4488: }

4490: /*@C
4491:   MatSeqAIJRestoreArrayRead - restore the read-only access array obtained from `MatSeqAIJGetArrayRead()`

4493:   Not Collective; No Fortran Support

4495:   Input Parameter:
4496: . A - a `MATSEQAIJ` matrix

4498:   Output Parameter:
4499: . array - pointer to the data

4501:   Level: intermediate

4503: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()`
4504: @*/
4505: PetscErrorCode MatSeqAIJRestoreArrayRead(Mat A, const PetscScalar *array[])
4506: {
4507:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4509:   PetscFunctionBegin;
4510:   if (aij->ops->restorearrayread) {
4511:     PetscCall((*aij->ops->restorearrayread)(A, array));
4512:   } else {
4513:     *array = NULL;
4514:   }
4515:   PetscFunctionReturn(PETSC_SUCCESS);
4516: }

4518: /*@C
4519:   MatSeqAIJGetArrayWrite - gives write-only access to the array where the data for a `MATSEQAIJ` matrix is stored

4521:   Not Collective; No Fortran Support

4523:   Input Parameter:
4524: . A - a `MATSEQAIJ` matrix

4526:   Output Parameter:
4527: . array - pointer to the data

4529:   Level: intermediate

4531: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayRead()`
4532: @*/
4533: PetscErrorCode MatSeqAIJGetArrayWrite(Mat A, PetscScalar *array[])
4534: {
4535:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4537:   PetscFunctionBegin;
4538:   if (aij->ops->getarraywrite) {
4539:     PetscCall((*aij->ops->getarraywrite)(A, array));
4540:   } else {
4541:     *array = aij->a;
4542:   }
4543:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
4544:   PetscCall(PetscObjectStateIncrease((PetscObject)A));
4545:   PetscFunctionReturn(PETSC_SUCCESS);
4546: }

4548: /*@C
4549:   MatSeqAIJRestoreArrayWrite - restore the read-only access array obtained from MatSeqAIJGetArrayRead

4551:   Not Collective; No Fortran Support

4553:   Input Parameter:
4554: . A - a MATSEQAIJ matrix

4556:   Output Parameter:
4557: . array - pointer to the data

4559:   Level: intermediate

4561: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()`
4562: @*/
4563: PetscErrorCode MatSeqAIJRestoreArrayWrite(Mat A, PetscScalar *array[])
4564: {
4565:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4567:   PetscFunctionBegin;
4568:   if (aij->ops->restorearraywrite) {
4569:     PetscCall((*aij->ops->restorearraywrite)(A, array));
4570:   } else {
4571:     *array = NULL;
4572:   }
4573:   PetscFunctionReturn(PETSC_SUCCESS);
4574: }

4576: /*@C
4577:   MatSeqAIJGetCSRAndMemType - Get the CSR arrays and the memory type of the `MATSEQAIJ` matrix

4579:   Not Collective; No Fortran Support

4581:   Input Parameter:
4582: . mat - a matrix of type `MATSEQAIJ` or its subclasses

4584:   Output Parameters:
4585: + i     - row map array of the matrix
4586: . j     - column index array of the matrix
4587: . a     - data array of the matrix
4588: - mtype - memory type of the arrays

4590:   Level: developer

4592:   Notes:
4593:   Any of the output parameters can be `NULL`, in which case the corresponding value is not returned.
4594:   If mat is a device matrix, the arrays are on the device. Otherwise, they are on the host.

4596:   One can call this routine on a preallocated but not assembled matrix to just get the memory of the CSR underneath the matrix.
4597:   If the matrix is assembled, the data array `a` is guaranteed to have the latest values of the matrix.

4599: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()`
4600: @*/
4601: PetscErrorCode MatSeqAIJGetCSRAndMemType(Mat mat, const PetscInt *i[], const PetscInt *j[], PetscScalar *a[], PetscMemType *mtype)
4602: {
4603:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;

4605:   PetscFunctionBegin;
4606:   PetscCheck(mat->preallocated, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "matrix is not preallocated");
4607:   if (aij->ops->getcsrandmemtype) {
4608:     PetscCall((*aij->ops->getcsrandmemtype)(mat, i, j, a, mtype));
4609:   } else {
4610:     if (i) *i = aij->i;
4611:     if (j) *j = aij->j;
4612:     if (a) *a = aij->a;
4613:     if (mtype) *mtype = PETSC_MEMTYPE_HOST;
4614:   }
4615:   PetscFunctionReturn(PETSC_SUCCESS);
4616: }

4618: /*@
4619:   MatSeqAIJGetMaxRowNonzeros - returns the maximum number of nonzeros in any row

4621:   Not Collective

4623:   Input Parameter:
4624: . A - a `MATSEQAIJ` matrix

4626:   Output Parameter:
4627: . nz - the maximum number of nonzeros in any row

4629:   Level: intermediate

4631: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRestoreArray()`, `MatSeqAIJGetArrayF90()`
4632: @*/
4633: PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat A, PetscInt *nz)
4634: {
4635:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4637:   PetscFunctionBegin;
4638:   *nz = aij->rmax;
4639:   PetscFunctionReturn(PETSC_SUCCESS);
4640: }

4642: static PetscErrorCode MatCOOStructDestroy_SeqAIJ(void **data)
4643: {
4644:   MatCOOStruct_SeqAIJ *coo = (MatCOOStruct_SeqAIJ *)*data;

4646:   PetscFunctionBegin;
4647:   PetscCall(PetscFree(coo->perm));
4648:   PetscCall(PetscFree(coo->jmap));
4649:   PetscCall(PetscFree(coo));
4650:   PetscFunctionReturn(PETSC_SUCCESS);
4651: }

4653: PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat mat, PetscCount coo_n, PetscInt coo_i[], PetscInt coo_j[])
4654: {
4655:   MPI_Comm             comm;
4656:   PetscInt            *i, *j;
4657:   PetscInt             M, N, row, iprev;
4658:   PetscCount           k, p, q, nneg, nnz, start, end; /* Index the coo array, so use PetscCount as their type */
4659:   PetscInt            *Ai;                             /* Change to PetscCount once we use it for row pointers */
4660:   PetscInt            *Aj;
4661:   PetscScalar         *Aa;
4662:   Mat_SeqAIJ          *seqaij = (Mat_SeqAIJ *)mat->data;
4663:   MatType              rtype;
4664:   PetscCount          *perm, *jmap;
4665:   MatCOOStruct_SeqAIJ *coo;
4666:   PetscBool            isorted;
4667:   PetscBool            hypre;
4668:   const char          *name;

4670:   PetscFunctionBegin;
4671:   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
4672:   PetscCall(MatGetSize(mat, &M, &N));
4673:   i = coo_i;
4674:   j = coo_j;
4675:   PetscCall(PetscMalloc1(coo_n, &perm));

4677:   /* Ignore entries with negative row or col indices; at the same time, check if i[] is already sorted (e.g., MatConvert_AlJ_HYPRE results in this case) */
4678:   isorted = PETSC_TRUE;
4679:   iprev   = PETSC_INT_MIN;
4680:   for (k = 0; k < coo_n; k++) {
4681:     if (j[k] < 0) i[k] = -1;
4682:     if (isorted) {
4683:       if (i[k] < iprev) isorted = PETSC_FALSE;
4684:       else iprev = i[k];
4685:     }
4686:     perm[k] = k;
4687:   }

4689:   /* Sort by row if not already */
4690:   if (!isorted) PetscCall(PetscSortIntWithIntCountArrayPair(coo_n, i, j, perm));

4692:   /* Advance k to the first row with a non-negative index */
4693:   for (k = 0; k < coo_n; k++)
4694:     if (i[k] >= 0) break;
4695:   nneg = k;
4696:   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 */
4697:   nnz = 0;                                          /* Total number of unique nonzeros to be counted */
4698:   jmap++;                                           /* Inc jmap by 1 for convenience */

4700:   PetscCall(PetscShmgetAllocateArray(M + 1, sizeof(PetscInt), (void **)&Ai)); /* CSR of A */
4701:   PetscCall(PetscArrayzero(Ai, M + 1));
4702:   PetscCall(PetscShmgetAllocateArray(coo_n - nneg, sizeof(PetscInt), (void **)&Aj)); /* We have at most coo_n-nneg unique nonzeros */

4704:   PetscCall(PetscObjectGetName((PetscObject)mat, &name));
4705:   PetscCall(PetscStrcmp("_internal_COO_mat_for_hypre", name, &hypre));

4707:   /* In each row, sort by column, then unique column indices to get row length */
4708:   Ai++;  /* Inc by 1 for convenience */
4709:   q = 0; /* q-th unique nonzero, with q starting from 0 */
4710:   while (k < coo_n) {
4711:     PetscBool strictly_sorted; // this row is strictly sorted?
4712:     PetscInt  jprev;

4714:     /* get [start,end) indices for this row; also check if cols in this row are strictly sorted */
4715:     row             = i[k];
4716:     start           = k;
4717:     jprev           = PETSC_INT_MIN;
4718:     strictly_sorted = PETSC_TRUE;
4719:     while (k < coo_n && i[k] == row) {
4720:       if (strictly_sorted) {
4721:         if (j[k] <= jprev) strictly_sorted = PETSC_FALSE;
4722:         else jprev = j[k];
4723:       }
4724:       k++;
4725:     }
4726:     end = k;

4728:     /* hack for HYPRE: swap min column to diag so that diagonal values will go first */
4729:     if (hypre) {
4730:       PetscInt  minj    = PETSC_INT_MAX;
4731:       PetscBool hasdiag = PETSC_FALSE;

4733:       if (strictly_sorted) { // fast path to swap the first and the diag
4734:         PetscCount tmp;
4735:         for (p = start; p < end; p++) {
4736:           if (j[p] == row && p != start) {
4737:             j[p]        = j[start];
4738:             j[start]    = row;
4739:             tmp         = perm[start];
4740:             perm[start] = perm[p];
4741:             perm[p]     = tmp;
4742:             break;
4743:           }
4744:         }
4745:       } else {
4746:         for (p = start; p < end; p++) {
4747:           hasdiag = (PetscBool)(hasdiag || (j[p] == row));
4748:           minj    = PetscMin(minj, j[p]);
4749:         }

4751:         if (hasdiag) {
4752:           for (p = start; p < end; p++) {
4753:             if (j[p] == minj) j[p] = row;
4754:             else if (j[p] == row) j[p] = minj;
4755:           }
4756:         }
4757:       }
4758:     }
4759:     // sort by columns in a row
4760:     if (!strictly_sorted) PetscCall(PetscSortIntWithCountArray(end - start, j + start, perm + start));

4762:     if (strictly_sorted) { // fast path to set Aj[], jmap[], Ai[], nnz, q
4763:       for (p = start; p < end; p++, q++) {
4764:         Aj[q]   = j[p];
4765:         jmap[q] = 1;
4766:       }
4767:       PetscCall(PetscIntCast(end - start, Ai + row));
4768:       nnz += Ai[row]; // q is already advanced
4769:     } else {
4770:       /* Find number of unique col entries in this row */
4771:       Aj[q]   = j[start]; /* Log the first nonzero in this row */
4772:       jmap[q] = 1;        /* Number of repeats of this nonzero entry */
4773:       Ai[row] = 1;
4774:       nnz++;

4776:       for (p = start + 1; p < end; p++) { /* Scan remaining nonzero in this row */
4777:         if (j[p] != j[p - 1]) {           /* Meet a new nonzero */
4778:           q++;
4779:           jmap[q] = 1;
4780:           Aj[q]   = j[p];
4781:           Ai[row]++;
4782:           nnz++;
4783:         } else {
4784:           jmap[q]++;
4785:         }
4786:       }
4787:       q++; /* Move to next row and thus next unique nonzero */
4788:     }
4789:   }

4791:   Ai--; /* Back to the beginning of Ai[] */
4792:   for (k = 0; k < M; k++) Ai[k + 1] += Ai[k];
4793:   jmap--; // Back to the beginning of jmap[]
4794:   jmap[0] = 0;
4795:   for (k = 0; k < nnz; k++) jmap[k + 1] += jmap[k];

4797:   if (nnz < coo_n - nneg) { /* Reallocate with actual number of unique nonzeros */
4798:     PetscCount *jmap_new;
4799:     PetscInt   *Aj_new;

4801:     PetscCall(PetscMalloc1(nnz + 1, &jmap_new));
4802:     PetscCall(PetscArraycpy(jmap_new, jmap, nnz + 1));
4803:     PetscCall(PetscFree(jmap));
4804:     jmap = jmap_new;

4806:     PetscCall(PetscShmgetAllocateArray(nnz, sizeof(PetscInt), (void **)&Aj_new));
4807:     PetscCall(PetscArraycpy(Aj_new, Aj, nnz));
4808:     PetscCall(PetscShmgetDeallocateArray((void **)&Aj));
4809:     Aj = Aj_new;
4810:   }

4812:   if (nneg) { /* Discard heading entries with negative indices in perm[], as we'll access it from index 0 in MatSetValuesCOO */
4813:     PetscCount *perm_new;

4815:     PetscCall(PetscMalloc1(coo_n - nneg, &perm_new));
4816:     PetscCall(PetscArraycpy(perm_new, perm + nneg, coo_n - nneg));
4817:     PetscCall(PetscFree(perm));
4818:     perm = perm_new;
4819:   }

4821:   PetscCall(MatGetRootType_Private(mat, &rtype));
4822:   PetscCall(PetscShmgetAllocateArray(nnz, sizeof(PetscScalar), (void **)&Aa));
4823:   PetscCall(PetscArrayzero(Aa, nnz));
4824:   PetscCall(MatSetSeqAIJWithArrays_private(PETSC_COMM_SELF, M, N, Ai, Aj, Aa, rtype, mat));

4826:   seqaij->free_a = seqaij->free_ij = PETSC_TRUE; /* Let newmat own Ai, Aj and Aa */

4828:   // Put the COO struct in a container and then attach that to the matrix
4829:   PetscCall(PetscMalloc1(1, &coo));
4830:   PetscCall(PetscIntCast(nnz, &coo->nz));
4831:   coo->n    = coo_n;
4832:   coo->Atot = coo_n - nneg; // Annz is seqaij->nz, so no need to record that again
4833:   coo->jmap = jmap;         // of length nnz+1
4834:   coo->perm = perm;
4835:   PetscCall(PetscObjectContainerCompose((PetscObject)mat, "__PETSc_MatCOOStruct_Host", coo, MatCOOStructDestroy_SeqAIJ));
4836:   PetscFunctionReturn(PETSC_SUCCESS);
4837: }

4839: static PetscErrorCode MatSetValuesCOO_SeqAIJ(Mat A, const PetscScalar v[], InsertMode imode)
4840: {
4841:   Mat_SeqAIJ          *aseq = (Mat_SeqAIJ *)A->data;
4842:   PetscCount           i, j, Annz = aseq->nz;
4843:   PetscCount          *perm, *jmap;
4844:   PetscScalar         *Aa;
4845:   PetscContainer       container;
4846:   MatCOOStruct_SeqAIJ *coo;

4848:   PetscFunctionBegin;
4849:   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_MatCOOStruct_Host", (PetscObject *)&container));
4850:   PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not found MatCOOStruct on this matrix");
4851:   PetscCall(PetscContainerGetPointer(container, (void **)&coo));
4852:   perm = coo->perm;
4853:   jmap = coo->jmap;
4854:   PetscCall(MatSeqAIJGetArray(A, &Aa));
4855:   for (i = 0; i < Annz; i++) {
4856:     PetscScalar sum = 0.0;
4857:     for (j = jmap[i]; j < jmap[i + 1]; j++) sum += v[perm[j]];
4858:     Aa[i] = (imode == INSERT_VALUES ? 0.0 : Aa[i]) + sum;
4859:   }
4860:   PetscCall(MatSeqAIJRestoreArray(A, &Aa));
4861:   PetscFunctionReturn(PETSC_SUCCESS);
4862: }

4864: #if defined(PETSC_HAVE_CUDA)
4865: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCUSPARSE(Mat, MatType, MatReuse, Mat *);
4866: #endif
4867: #if defined(PETSC_HAVE_HIP)
4868: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJHIPSPARSE(Mat, MatType, MatReuse, Mat *);
4869: #endif
4870: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
4871: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJKokkos(Mat, MatType, MatReuse, Mat *);
4872: #endif

4874: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat B)
4875: {
4876:   Mat_SeqAIJ *b;
4877:   PetscMPIInt size;

4879:   PetscFunctionBegin;
4880:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
4881:   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Comm must be of size 1");

4883:   PetscCall(PetscNew(&b));

4885:   B->data   = (void *)b;
4886:   B->ops[0] = MatOps_Values;
4887:   if (B->sortedfull) B->ops->setvalues = MatSetValues_SeqAIJ_SortedFull;

4889:   b->row                = NULL;
4890:   b->col                = NULL;
4891:   b->icol               = NULL;
4892:   b->reallocs           = 0;
4893:   b->ignorezeroentries  = PETSC_FALSE;
4894:   b->roworiented        = PETSC_TRUE;
4895:   b->nonew              = 0;
4896:   b->diag               = NULL;
4897:   b->solve_work         = NULL;
4898:   B->spptr              = NULL;
4899:   b->saved_values       = NULL;
4900:   b->idiag              = NULL;
4901:   b->mdiag              = NULL;
4902:   b->ssor_work          = NULL;
4903:   b->omega              = 1.0;
4904:   b->fshift             = 0.0;
4905:   b->idiagvalid         = PETSC_FALSE;
4906:   b->ibdiagvalid        = PETSC_FALSE;
4907:   b->keepnonzeropattern = PETSC_FALSE;

4909:   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ));
4910: #if defined(PETSC_HAVE_MATLAB)
4911:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "PetscMatlabEnginePut_C", MatlabEnginePut_SeqAIJ));
4912:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "PetscMatlabEngineGet_C", MatlabEngineGet_SeqAIJ));
4913: #endif
4914:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJSetColumnIndices_C", MatSeqAIJSetColumnIndices_SeqAIJ));
4915:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_SeqAIJ));
4916:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_SeqAIJ));
4917:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqsbaij_C", MatConvert_SeqAIJ_SeqSBAIJ));
4918:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqbaij_C", MatConvert_SeqAIJ_SeqBAIJ));
4919:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijperm_C", MatConvert_SeqAIJ_SeqAIJPERM));
4920:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijsell_C", MatConvert_SeqAIJ_SeqAIJSELL));
4921: #if defined(PETSC_HAVE_MKL_SPARSE)
4922:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijmkl_C", MatConvert_SeqAIJ_SeqAIJMKL));
4923: #endif
4924: #if defined(PETSC_HAVE_CUDA)
4925:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijcusparse_C", MatConvert_SeqAIJ_SeqAIJCUSPARSE));
4926:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaijcusparse_seqaij_C", MatProductSetFromOptions_SeqAIJ));
4927:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqaijcusparse_C", MatProductSetFromOptions_SeqAIJ));
4928: #endif
4929: #if defined(PETSC_HAVE_HIP)
4930:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijhipsparse_C", MatConvert_SeqAIJ_SeqAIJHIPSPARSE));
4931:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaijhipsparse_seqaij_C", MatProductSetFromOptions_SeqAIJ));
4932:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqaijhipsparse_C", MatProductSetFromOptions_SeqAIJ));
4933: #endif
4934: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
4935:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijkokkos_C", MatConvert_SeqAIJ_SeqAIJKokkos));
4936: #endif
4937:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqaijcrl_C", MatConvert_SeqAIJ_SeqAIJCRL));
4938: #if defined(PETSC_HAVE_ELEMENTAL)
4939:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_elemental_C", MatConvert_SeqAIJ_Elemental));
4940: #endif
4941: #if defined(PETSC_HAVE_SCALAPACK)
4942:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_scalapack_C", MatConvert_AIJ_ScaLAPACK));
4943: #endif
4944: #if defined(PETSC_HAVE_HYPRE)
4945:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_hypre_C", MatConvert_AIJ_HYPRE));
4946:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_transpose_seqaij_seqaij_C", MatProductSetFromOptions_Transpose_AIJ_AIJ));
4947: #endif
4948:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqdense_C", MatConvert_SeqAIJ_SeqDense));
4949:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_seqsell_C", MatConvert_SeqAIJ_SeqSELL));
4950:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaij_is_C", MatConvert_XAIJ_IS));
4951:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsTranspose_C", MatIsTranspose_SeqAIJ));
4952:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsHermitianTranspose_C", MatIsHermitianTranspose_SeqAIJ));
4953:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJSetPreallocation_C", MatSeqAIJSetPreallocation_SeqAIJ));
4954:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatResetPreallocation_C", MatResetPreallocation_SeqAIJ));
4955:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatResetHash_C", MatResetHash_SeqAIJ));
4956:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJSetPreallocationCSR_C", MatSeqAIJSetPreallocationCSR_SeqAIJ));
4957:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatReorderForNonzeroDiagonal_C", MatReorderForNonzeroDiagonal_SeqAIJ));
4958:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_is_seqaij_C", MatProductSetFromOptions_IS_XAIJ));
4959:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqdense_seqaij_C", MatProductSetFromOptions_SeqDense_SeqAIJ));
4960:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqaij_C", MatProductSetFromOptions_SeqAIJ));
4961:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqAIJKron_C", MatSeqAIJKron_SeqAIJ));
4962:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSetPreallocationCOO_C", MatSetPreallocationCOO_SeqAIJ));
4963:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSetValuesCOO_C", MatSetValuesCOO_SeqAIJ));
4964:   PetscCall(MatCreate_SeqAIJ_Inode(B));
4965:   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ));
4966:   PetscCall(MatSeqAIJSetTypeFromOptions(B)); /* this allows changing the matrix subtype to say MATSEQAIJPERM */
4967:   PetscFunctionReturn(PETSC_SUCCESS);
4968: }

4970: /*
4971:     Given a matrix generated with MatGetFactor() duplicates all the information in A into C
4972: */
4973: PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat C, Mat A, MatDuplicateOption cpvalues, PetscBool mallocmatspace)
4974: {
4975:   Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data, *a = (Mat_SeqAIJ *)A->data;
4976:   PetscInt    m = A->rmap->n, i;

4978:   PetscFunctionBegin;
4979:   PetscCheck(A->assembled || cpvalues == MAT_DO_NOT_COPY_VALUES, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot duplicate unassembled matrix");

4981:   C->factortype    = A->factortype;
4982:   c->row           = NULL;
4983:   c->col           = NULL;
4984:   c->icol          = NULL;
4985:   c->reallocs      = 0;
4986:   c->diagonaldense = a->diagonaldense;

4988:   C->assembled = A->assembled;

4990:   if (A->preallocated) {
4991:     PetscCall(PetscLayoutReference(A->rmap, &C->rmap));
4992:     PetscCall(PetscLayoutReference(A->cmap, &C->cmap));

4994:     if (!A->hash_active) {
4995:       PetscCall(PetscMalloc1(m, &c->imax));
4996:       PetscCall(PetscMemcpy(c->imax, a->imax, m * sizeof(PetscInt)));
4997:       PetscCall(PetscMalloc1(m, &c->ilen));
4998:       PetscCall(PetscMemcpy(c->ilen, a->ilen, m * sizeof(PetscInt)));

5000:       /* allocate the matrix space */
5001:       if (mallocmatspace) {
5002:         PetscCall(PetscShmgetAllocateArray(a->i[m], sizeof(PetscScalar), (void **)&c->a));
5003:         PetscCall(PetscShmgetAllocateArray(a->i[m], sizeof(PetscInt), (void **)&c->j));
5004:         PetscCall(PetscShmgetAllocateArray(m + 1, sizeof(PetscInt), (void **)&c->i));
5005:         PetscCall(PetscArraycpy(c->i, a->i, m + 1));
5006:         c->free_a  = PETSC_TRUE;
5007:         c->free_ij = PETSC_TRUE;
5008:         if (m > 0) {
5009:           PetscCall(PetscArraycpy(c->j, a->j, a->i[m]));
5010:           if (cpvalues == MAT_COPY_VALUES) {
5011:             const PetscScalar *aa;

5013:             PetscCall(MatSeqAIJGetArrayRead(A, &aa));
5014:             PetscCall(PetscArraycpy(c->a, aa, a->i[m]));
5015:             PetscCall(MatSeqAIJGetArrayRead(A, &aa));
5016:           } else {
5017:             PetscCall(PetscArrayzero(c->a, a->i[m]));
5018:           }
5019:         }
5020:       }
5021:       C->preallocated = PETSC_TRUE;
5022:     } else {
5023:       PetscCheck(mallocmatspace, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Cannot malloc matrix memory from a non-preallocated matrix");
5024:       PetscCall(MatSetUp(C));
5025:     }

5027:     c->ignorezeroentries = a->ignorezeroentries;
5028:     c->roworiented       = a->roworiented;
5029:     c->nonew             = a->nonew;
5030:     if (a->diag) {
5031:       PetscCall(PetscMalloc1(m + 1, &c->diag));
5032:       PetscCall(PetscMemcpy(c->diag, a->diag, m * sizeof(PetscInt)));
5033:     } else c->diag = NULL;

5035:     c->solve_work         = NULL;
5036:     c->saved_values       = NULL;
5037:     c->idiag              = NULL;
5038:     c->ssor_work          = NULL;
5039:     c->keepnonzeropattern = a->keepnonzeropattern;

5041:     c->rmax  = a->rmax;
5042:     c->nz    = a->nz;
5043:     c->maxnz = a->nz; /* Since we allocate exactly the right amount */

5045:     c->compressedrow.use   = a->compressedrow.use;
5046:     c->compressedrow.nrows = a->compressedrow.nrows;
5047:     if (a->compressedrow.use) {
5048:       i = a->compressedrow.nrows;
5049:       PetscCall(PetscMalloc2(i + 1, &c->compressedrow.i, i, &c->compressedrow.rindex));
5050:       PetscCall(PetscArraycpy(c->compressedrow.i, a->compressedrow.i, i + 1));
5051:       PetscCall(PetscArraycpy(c->compressedrow.rindex, a->compressedrow.rindex, i));
5052:     } else {
5053:       c->compressedrow.use    = PETSC_FALSE;
5054:       c->compressedrow.i      = NULL;
5055:       c->compressedrow.rindex = NULL;
5056:     }
5057:     c->nonzerorowcnt = a->nonzerorowcnt;
5058:     C->nonzerostate  = A->nonzerostate;

5060:     PetscCall(MatDuplicate_SeqAIJ_Inode(A, cpvalues, &C));
5061:   }
5062:   PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist, &((PetscObject)C)->qlist));
5063:   PetscFunctionReturn(PETSC_SUCCESS);
5064: }

5066: PetscErrorCode MatDuplicate_SeqAIJ(Mat A, MatDuplicateOption cpvalues, Mat *B)
5067: {
5068:   PetscFunctionBegin;
5069:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
5070:   PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n));
5071:   if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) PetscCall(MatSetBlockSizesFromMats(*B, A, A));
5072:   PetscCall(MatSetType(*B, ((PetscObject)A)->type_name));
5073:   PetscCall(MatDuplicateNoCreate_SeqAIJ(*B, A, cpvalues, PETSC_TRUE));
5074:   PetscFunctionReturn(PETSC_SUCCESS);
5075: }

5077: PetscErrorCode MatLoad_SeqAIJ(Mat newMat, PetscViewer viewer)
5078: {
5079:   PetscBool isbinary, ishdf5;

5081:   PetscFunctionBegin;
5084:   /* force binary viewer to load .info file if it has not yet done so */
5085:   PetscCall(PetscViewerSetUp(viewer));
5086:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
5087:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
5088:   if (isbinary) {
5089:     PetscCall(MatLoad_SeqAIJ_Binary(newMat, viewer));
5090:   } else if (ishdf5) {
5091: #if defined(PETSC_HAVE_HDF5)
5092:     PetscCall(MatLoad_AIJ_HDF5(newMat, viewer));
5093: #else
5094:     SETERRQ(PetscObjectComm((PetscObject)newMat), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
5095: #endif
5096:   } else {
5097:     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);
5098:   }
5099:   PetscFunctionReturn(PETSC_SUCCESS);
5100: }

5102: PetscErrorCode MatLoad_SeqAIJ_Binary(Mat mat, PetscViewer viewer)
5103: {
5104:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)mat->data;
5105:   PetscInt    header[4], *rowlens, M, N, nz, sum, rows, cols, i;

5107:   PetscFunctionBegin;
5108:   PetscCall(PetscViewerSetUp(viewer));

5110:   /* read in matrix header */
5111:   PetscCall(PetscViewerBinaryRead(viewer, header, 4, NULL, PETSC_INT));
5112:   PetscCheck(header[0] == MAT_FILE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a matrix object in file");
5113:   M  = header[1];
5114:   N  = header[2];
5115:   nz = header[3];
5116:   PetscCheck(M >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Matrix row size (%" PetscInt_FMT ") in file is negative", M);
5117:   PetscCheck(N >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Matrix column size (%" PetscInt_FMT ") in file is negative", N);
5118:   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Matrix stored in special format on disk, cannot load as SeqAIJ");

5120:   /* set block sizes from the viewer's .info file */
5121:   PetscCall(MatLoad_Binary_BlockSizes(mat, viewer));
5122:   /* set local and global sizes if not set already */
5123:   if (mat->rmap->n < 0) mat->rmap->n = M;
5124:   if (mat->cmap->n < 0) mat->cmap->n = N;
5125:   if (mat->rmap->N < 0) mat->rmap->N = M;
5126:   if (mat->cmap->N < 0) mat->cmap->N = N;
5127:   PetscCall(PetscLayoutSetUp(mat->rmap));
5128:   PetscCall(PetscLayoutSetUp(mat->cmap));

5130:   /* check if the matrix sizes are correct */
5131:   PetscCall(MatGetSize(mat, &rows, &cols));
5132:   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);

5134:   /* read in row lengths */
5135:   PetscCall(PetscMalloc1(M, &rowlens));
5136:   PetscCall(PetscViewerBinaryRead(viewer, rowlens, M, NULL, PETSC_INT));
5137:   /* check if sum(rowlens) is same as nz */
5138:   sum = 0;
5139:   for (i = 0; i < M; i++) sum += rowlens[i];
5140:   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);
5141:   /* preallocate and check sizes */
5142:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(mat, 0, rowlens));
5143:   PetscCall(MatGetSize(mat, &rows, &cols));
5144:   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);
5145:   /* store row lengths */
5146:   PetscCall(PetscArraycpy(a->ilen, rowlens, M));
5147:   PetscCall(PetscFree(rowlens));

5149:   /* fill in "i" row pointers */
5150:   a->i[0] = 0;
5151:   for (i = 0; i < M; i++) a->i[i + 1] = a->i[i] + a->ilen[i];
5152:   /* read in "j" column indices */
5153:   PetscCall(PetscViewerBinaryRead(viewer, a->j, nz, NULL, PETSC_INT));
5154:   /* read in "a" nonzero values */
5155:   PetscCall(PetscViewerBinaryRead(viewer, a->a, nz, NULL, PETSC_SCALAR));

5157:   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
5158:   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
5159:   PetscFunctionReturn(PETSC_SUCCESS);
5160: }

5162: PetscErrorCode MatEqual_SeqAIJ(Mat A, Mat B, PetscBool *flg)
5163: {
5164:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
5165:   const PetscScalar *aa, *ba;
5166: #if defined(PETSC_USE_COMPLEX)
5167:   PetscInt k;
5168: #endif

5170:   PetscFunctionBegin;
5171:   /* If the  matrix dimensions are not equal,or no of nonzeros */
5172:   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) || (a->nz != b->nz)) {
5173:     *flg = PETSC_FALSE;
5174:     PetscFunctionReturn(PETSC_SUCCESS);
5175:   }

5177:   /* if the a->i are the same */
5178:   PetscCall(PetscArraycmp(a->i, b->i, A->rmap->n + 1, flg));
5179:   if (!*flg) PetscFunctionReturn(PETSC_SUCCESS);

5181:   /* if a->j are the same */
5182:   PetscCall(PetscArraycmp(a->j, b->j, a->nz, flg));
5183:   if (!*flg) PetscFunctionReturn(PETSC_SUCCESS);

5185:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
5186:   PetscCall(MatSeqAIJGetArrayRead(B, &ba));
5187:   /* if a->a are the same */
5188: #if defined(PETSC_USE_COMPLEX)
5189:   for (k = 0; k < a->nz; k++) {
5190:     if (PetscRealPart(aa[k]) != PetscRealPart(ba[k]) || PetscImaginaryPart(aa[k]) != PetscImaginaryPart(ba[k])) {
5191:       *flg = PETSC_FALSE;
5192:       PetscFunctionReturn(PETSC_SUCCESS);
5193:     }
5194:   }
5195: #else
5196:   PetscCall(PetscArraycmp(aa, ba, a->nz, flg));
5197: #endif
5198:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
5199:   PetscCall(MatSeqAIJRestoreArrayRead(B, &ba));
5200:   PetscFunctionReturn(PETSC_SUCCESS);
5201: }

5203: /*@
5204:   MatCreateSeqAIJWithArrays - Creates an sequential `MATSEQAIJ` matrix using matrix elements (in CSR format)
5205:   provided by the user.

5207:   Collective

5209:   Input Parameters:
5210: + comm - must be an MPI communicator of size 1
5211: . m    - number of rows
5212: . n    - number of columns
5213: . i    - row indices; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
5214: . j    - column indices
5215: - a    - matrix values

5217:   Output Parameter:
5218: . mat - the matrix

5220:   Level: intermediate

5222:   Notes:
5223:   The `i`, `j`, and `a` arrays are not copied by this routine, the user must free these arrays
5224:   once the matrix is destroyed and not before

5226:   You cannot set new nonzero locations into this matrix, that will generate an error.

5228:   The `i` and `j` indices are 0 based

5230:   The format which is used for the sparse matrix input, is equivalent to a
5231:   row-major ordering.. i.e for the following matrix, the input data expected is
5232:   as shown
5233: .vb
5234:         1 0 0
5235:         2 0 3
5236:         4 5 6

5238:         i =  {0,1,3,6}  [size = nrow+1  = 3+1]
5239:         j =  {0,0,2,0,1,2}  [size = 6]; values must be sorted for each row
5240:         v =  {1,2,3,4,5,6}  [size = 6]
5241: .ve

5243: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MatCreateMPIAIJWithArrays()`, `MatMPIAIJSetPreallocationCSR()`
5244: @*/
5245: PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt i[], PetscInt j[], PetscScalar a[], Mat *mat)
5246: {
5247:   PetscInt    ii;
5248:   Mat_SeqAIJ *aij;
5249:   PetscInt    jj;

5251:   PetscFunctionBegin;
5252:   PetscCheck(m <= 0 || i[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "i (row indices) must start with 0");
5253:   PetscCall(MatCreate(comm, mat));
5254:   PetscCall(MatSetSizes(*mat, m, n, m, n));
5255:   /* PetscCall(MatSetBlockSizes(*mat,,)); */
5256:   PetscCall(MatSetType(*mat, MATSEQAIJ));
5257:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*mat, MAT_SKIP_ALLOCATION, NULL));
5258:   aij = (Mat_SeqAIJ *)(*mat)->data;
5259:   PetscCall(PetscMalloc1(m, &aij->imax));
5260:   PetscCall(PetscMalloc1(m, &aij->ilen));

5262:   aij->i       = i;
5263:   aij->j       = j;
5264:   aij->a       = a;
5265:   aij->nonew   = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
5266:   aij->free_a  = PETSC_FALSE;
5267:   aij->free_ij = PETSC_FALSE;

5269:   for (ii = 0, aij->nonzerorowcnt = 0, aij->rmax = 0; ii < m; ii++) {
5270:     aij->ilen[ii] = aij->imax[ii] = i[ii + 1] - i[ii];
5271:     if (PetscDefined(USE_DEBUG)) {
5272:       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]);
5273:       for (jj = i[ii] + 1; jj < i[ii + 1]; jj++) {
5274:         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);
5275:         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);
5276:       }
5277:     }
5278:   }
5279:   if (PetscDefined(USE_DEBUG)) {
5280:     for (ii = 0; ii < aij->i[m]; ii++) {
5281:       PetscCheck(j[ii] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative column index at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]);
5282:       PetscCheck(j[ii] <= n - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column index to large at location = %" PetscInt_FMT " index = %" PetscInt_FMT " last column = %" PetscInt_FMT, ii, j[ii], n - 1);
5283:     }
5284:   }

5286:   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
5287:   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
5288:   PetscFunctionReturn(PETSC_SUCCESS);
5289: }

5291: /*@
5292:   MatCreateSeqAIJFromTriple - Creates an sequential `MATSEQAIJ` matrix using matrix elements (in COO format)
5293:   provided by the user.

5295:   Collective

5297:   Input Parameters:
5298: + comm - must be an MPI communicator of size 1
5299: . m    - number of rows
5300: . n    - number of columns
5301: . i    - row indices
5302: . j    - column indices
5303: . a    - matrix values
5304: . nz   - number of nonzeros
5305: - idx  - if the `i` and `j` indices start with 1 use `PETSC_TRUE` otherwise use `PETSC_FALSE`

5307:   Output Parameter:
5308: . mat - the matrix

5310:   Level: intermediate

5312:   Example:
5313:   For the following matrix, the input data expected is as shown (using 0 based indexing)
5314: .vb
5315:         1 0 0
5316:         2 0 3
5317:         4 5 6

5319:         i =  {0,1,1,2,2,2}
5320:         j =  {0,0,2,0,1,2}
5321:         v =  {1,2,3,4,5,6}
5322: .ve

5324:   Note:
5325:   Instead of using this function, users should also consider `MatSetPreallocationCOO()` and `MatSetValuesCOO()`, which allow repeated or remote entries,
5326:   and are particularly useful in iterative applications.

5328: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MatCreateSeqAIJWithArrays()`, `MatMPIAIJSetPreallocationCSR()`, `MatSetValuesCOO()`, `MatSetPreallocationCOO()`
5329: @*/
5330: PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt i[], PetscInt j[], PetscScalar a[], Mat *mat, PetscInt nz, PetscBool idx)
5331: {
5332:   PetscInt ii, *nnz, one = 1, row, col;

5334:   PetscFunctionBegin;
5335:   PetscCall(PetscCalloc1(m, &nnz));
5336:   for (ii = 0; ii < nz; ii++) nnz[i[ii] - !!idx] += 1;
5337:   PetscCall(MatCreate(comm, mat));
5338:   PetscCall(MatSetSizes(*mat, m, n, m, n));
5339:   PetscCall(MatSetType(*mat, MATSEQAIJ));
5340:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*mat, 0, nnz));
5341:   for (ii = 0; ii < nz; ii++) {
5342:     if (idx) {
5343:       row = i[ii] - 1;
5344:       col = j[ii] - 1;
5345:     } else {
5346:       row = i[ii];
5347:       col = j[ii];
5348:     }
5349:     PetscCall(MatSetValues(*mat, one, &row, one, &col, &a[ii], ADD_VALUES));
5350:   }
5351:   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
5352:   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
5353:   PetscCall(PetscFree(nnz));
5354:   PetscFunctionReturn(PETSC_SUCCESS);
5355: }

5357: PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat A)
5358: {
5359:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

5361:   PetscFunctionBegin;
5362:   a->idiagvalid  = PETSC_FALSE;
5363:   a->ibdiagvalid = PETSC_FALSE;

5365:   PetscCall(MatSeqAIJInvalidateDiagonal_Inode(A));
5366:   PetscFunctionReturn(PETSC_SUCCESS);
5367: }

5369: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm comm, Mat inmat, PetscInt n, MatReuse scall, Mat *outmat)
5370: {
5371:   PetscFunctionBegin;
5372:   PetscCall(MatCreateMPIMatConcatenateSeqMat_MPIAIJ(comm, inmat, n, scall, outmat));
5373:   PetscFunctionReturn(PETSC_SUCCESS);
5374: }

5376: /*
5377:  Permute A into C's *local* index space using rowemb,colemb.
5378:  The embedding are supposed to be injections and the above implies that the range of rowemb is a subset
5379:  of [0,m), colemb is in [0,n).
5380:  If pattern == DIFFERENT_NONZERO_PATTERN, C is preallocated according to A.
5381:  */
5382: PetscErrorCode MatSetSeqMat_SeqAIJ(Mat C, IS rowemb, IS colemb, MatStructure pattern, Mat B)
5383: {
5384:   /* If making this function public, change the error returned in this function away from _PLIB. */
5385:   Mat_SeqAIJ     *Baij;
5386:   PetscBool       seqaij;
5387:   PetscInt        m, n, *nz, i, j, count;
5388:   PetscScalar     v;
5389:   const PetscInt *rowindices, *colindices;

5391:   PetscFunctionBegin;
5392:   if (!B) PetscFunctionReturn(PETSC_SUCCESS);
5393:   /* Check to make sure the target matrix (and embeddings) are compatible with C and each other. */
5394:   PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATSEQAIJ, &seqaij));
5395:   PetscCheck(seqaij, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is of wrong type");
5396:   if (rowemb) {
5397:     PetscCall(ISGetLocalSize(rowemb, &m));
5398:     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);
5399:   } else {
5400:     PetscCheck(C->rmap->n == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is row-incompatible with the target matrix");
5401:   }
5402:   if (colemb) {
5403:     PetscCall(ISGetLocalSize(colemb, &n));
5404:     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);
5405:   } else {
5406:     PetscCheck(C->cmap->n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is col-incompatible with the target matrix");
5407:   }

5409:   Baij = (Mat_SeqAIJ *)B->data;
5410:   if (pattern == DIFFERENT_NONZERO_PATTERN) {
5411:     PetscCall(PetscMalloc1(B->rmap->n, &nz));
5412:     for (i = 0; i < B->rmap->n; i++) nz[i] = Baij->i[i + 1] - Baij->i[i];
5413:     PetscCall(MatSeqAIJSetPreallocation(C, 0, nz));
5414:     PetscCall(PetscFree(nz));
5415:   }
5416:   if (pattern == SUBSET_NONZERO_PATTERN) PetscCall(MatZeroEntries(C));
5417:   count      = 0;
5418:   rowindices = NULL;
5419:   colindices = NULL;
5420:   if (rowemb) PetscCall(ISGetIndices(rowemb, &rowindices));
5421:   if (colemb) PetscCall(ISGetIndices(colemb, &colindices));
5422:   for (i = 0; i < B->rmap->n; i++) {
5423:     PetscInt row;
5424:     row = i;
5425:     if (rowindices) row = rowindices[i];
5426:     for (j = Baij->i[i]; j < Baij->i[i + 1]; j++) {
5427:       PetscInt col;
5428:       col = Baij->j[count];
5429:       if (colindices) col = colindices[col];
5430:       v = Baij->a[count];
5431:       PetscCall(MatSetValues(C, 1, &row, 1, &col, &v, INSERT_VALUES));
5432:       ++count;
5433:     }
5434:   }
5435:   /* FIXME: set C's nonzerostate correctly. */
5436:   /* Assembly for C is necessary. */
5437:   C->preallocated  = PETSC_TRUE;
5438:   C->assembled     = PETSC_TRUE;
5439:   C->was_assembled = PETSC_FALSE;
5440:   PetscFunctionReturn(PETSC_SUCCESS);
5441: }

5443: PetscErrorCode MatEliminateZeros_SeqAIJ(Mat A, PetscBool keep)
5444: {
5445:   Mat_SeqAIJ *a  = (Mat_SeqAIJ *)A->data;
5446:   MatScalar  *aa = a->a;
5447:   PetscInt    m = A->rmap->n, fshift = 0, fshift_prev = 0, i, k;
5448:   PetscInt   *ailen = a->ilen, *imax = a->imax, *ai = a->i, *aj = a->j, rmax = 0;

5450:   PetscFunctionBegin;
5451:   PetscCheck(A->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot eliminate zeros for unassembled matrix");
5452:   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
5453:   for (i = 1; i <= m; i++) {
5454:     /* move each nonzero entry back by the amount of zero slots (fshift) before it*/
5455:     for (k = ai[i - 1]; k < ai[i]; k++) {
5456:       if (aa[k] == 0 && (aj[k] != i - 1 || !keep)) fshift++;
5457:       else {
5458:         if (aa[k] == 0 && aj[k] == i - 1) PetscCall(PetscInfo(A, "Keep the diagonal zero at row %" PetscInt_FMT "\n", i - 1));
5459:         aa[k - fshift] = aa[k];
5460:         aj[k - fshift] = aj[k];
5461:       }
5462:     }
5463:     ai[i - 1] -= fshift_prev; // safe to update ai[i-1] now since it will not be used in the next iteration
5464:     fshift_prev = fshift;
5465:     /* reset ilen and imax for each row */
5466:     ailen[i - 1] = imax[i - 1] = ai[i] - fshift - ai[i - 1];
5467:     a->nonzerorowcnt += ((ai[i] - fshift - ai[i - 1]) > 0);
5468:     rmax = PetscMax(rmax, ailen[i - 1]);
5469:   }
5470:   if (fshift) {
5471:     if (m) {
5472:       ai[m] -= fshift;
5473:       a->nz = ai[m];
5474:     }
5475:     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));
5476:     A->nonzerostate++;
5477:     A->info.nz_unneeded += (PetscReal)fshift;
5478:     a->rmax = rmax;
5479:     if (a->inode.use && a->inode.checked) PetscCall(MatSeqAIJCheckInode(A));
5480:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5481:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5482:   }
5483:   PetscFunctionReturn(PETSC_SUCCESS);
5484: }

5486: PetscFunctionList MatSeqAIJList = NULL;

5488: /*@
5489:   MatSeqAIJSetType - Converts a `MATSEQAIJ` matrix to a subtype

5491:   Collective

5493:   Input Parameters:
5494: + mat    - the matrix object
5495: - matype - matrix type

5497:   Options Database Key:
5498: . -mat_seqaij_type  <method> - for example seqaijcrl

5500:   Level: intermediate

5502: .seealso: [](ch_matrices), `Mat`, `PCSetType()`, `VecSetType()`, `MatCreate()`, `MatType`
5503: @*/
5504: PetscErrorCode MatSeqAIJSetType(Mat mat, MatType matype)
5505: {
5506:   PetscBool sametype;
5507:   PetscErrorCode (*r)(Mat, MatType, MatReuse, Mat *);

5509:   PetscFunctionBegin;
5511:   PetscCall(PetscObjectTypeCompare((PetscObject)mat, matype, &sametype));
5512:   if (sametype) PetscFunctionReturn(PETSC_SUCCESS);

5514:   PetscCall(PetscFunctionListFind(MatSeqAIJList, matype, &r));
5515:   PetscCheck(r, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Mat type given: %s", matype);
5516:   PetscCall((*r)(mat, matype, MAT_INPLACE_MATRIX, &mat));
5517:   PetscFunctionReturn(PETSC_SUCCESS);
5518: }

5520: /*@C
5521:   MatSeqAIJRegister -  - Adds a new sub-matrix type for sequential `MATSEQAIJ` matrices

5523:   Not Collective, No Fortran Support

5525:   Input Parameters:
5526: + sname    - name of a new user-defined matrix type, for example `MATSEQAIJCRL`
5527: - function - routine to convert to subtype

5529:   Level: advanced

5531:   Notes:
5532:   `MatSeqAIJRegister()` may be called multiple times to add several user-defined solvers.

5534:   Then, your matrix can be chosen with the procedural interface at runtime via the option
5535: $     -mat_seqaij_type my_mat

5537: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRegisterAll()`
5538: @*/
5539: PetscErrorCode MatSeqAIJRegister(const char sname[], PetscErrorCode (*function)(Mat, MatType, MatReuse, Mat *))
5540: {
5541:   PetscFunctionBegin;
5542:   PetscCall(MatInitializePackage());
5543:   PetscCall(PetscFunctionListAdd(&MatSeqAIJList, sname, function));
5544:   PetscFunctionReturn(PETSC_SUCCESS);
5545: }

5547: PetscBool MatSeqAIJRegisterAllCalled = PETSC_FALSE;

5549: /*@C
5550:   MatSeqAIJRegisterAll - Registers all of the matrix subtypes of `MATSSEQAIJ`

5552:   Not Collective

5554:   Level: advanced

5556:   Note:
5557:   This registers the versions of `MATSEQAIJ` for GPUs

5559: .seealso: [](ch_matrices), `Mat`, `MatRegisterAll()`, `MatSeqAIJRegister()`
5560: @*/
5561: PetscErrorCode MatSeqAIJRegisterAll(void)
5562: {
5563:   PetscFunctionBegin;
5564:   if (MatSeqAIJRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
5565:   MatSeqAIJRegisterAllCalled = PETSC_TRUE;

5567:   PetscCall(MatSeqAIJRegister(MATSEQAIJCRL, MatConvert_SeqAIJ_SeqAIJCRL));
5568:   PetscCall(MatSeqAIJRegister(MATSEQAIJPERM, MatConvert_SeqAIJ_SeqAIJPERM));
5569:   PetscCall(MatSeqAIJRegister(MATSEQAIJSELL, MatConvert_SeqAIJ_SeqAIJSELL));
5570: #if defined(PETSC_HAVE_MKL_SPARSE)
5571:   PetscCall(MatSeqAIJRegister(MATSEQAIJMKL, MatConvert_SeqAIJ_SeqAIJMKL));
5572: #endif
5573: #if defined(PETSC_HAVE_CUDA)
5574:   PetscCall(MatSeqAIJRegister(MATSEQAIJCUSPARSE, MatConvert_SeqAIJ_SeqAIJCUSPARSE));
5575: #endif
5576: #if defined(PETSC_HAVE_HIP)
5577:   PetscCall(MatSeqAIJRegister(MATSEQAIJHIPSPARSE, MatConvert_SeqAIJ_SeqAIJHIPSPARSE));
5578: #endif
5579: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
5580:   PetscCall(MatSeqAIJRegister(MATSEQAIJKOKKOS, MatConvert_SeqAIJ_SeqAIJKokkos));
5581: #endif
5582: #if defined(PETSC_HAVE_VIENNACL) && defined(PETSC_HAVE_VIENNACL_NO_CUDA)
5583:   PetscCall(MatSeqAIJRegister(MATMPIAIJVIENNACL, MatConvert_SeqAIJ_SeqAIJViennaCL));
5584: #endif
5585:   PetscFunctionReturn(PETSC_SUCCESS);
5586: }

5588: /*
5589:     Special version for direct calls from Fortran
5590: */
5591: #if defined(PETSC_HAVE_FORTRAN_CAPS)
5592:   #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
5593: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
5594:   #define matsetvaluesseqaij_ matsetvaluesseqaij
5595: #endif

5597: /* Change these macros so can be used in void function */

5599: /* Change these macros so can be used in void function */
5600: /* Identical to PetscCallVoid, except it assigns to *_ierr */
5601: #undef PetscCall
5602: #define PetscCall(...) \
5603:   do { \
5604:     PetscErrorCode ierr_msv_mpiaij = __VA_ARGS__; \
5605:     if (PetscUnlikely(ierr_msv_mpiaij)) { \
5606:       *_ierr = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_msv_mpiaij, PETSC_ERROR_REPEAT, " "); \
5607:       return; \
5608:     } \
5609:   } while (0)

5611: #undef SETERRQ
5612: #define SETERRQ(comm, ierr, ...) \
5613:   do { \
5614:     *_ierr = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
5615:     return; \
5616:   } while (0)

5618: PETSC_EXTERN void matsetvaluesseqaij_(Mat *AA, PetscInt *mm, const PetscInt im[], PetscInt *nn, const PetscInt in[], const PetscScalar v[], InsertMode *isis, PetscErrorCode *_ierr)
5619: {
5620:   Mat         A = *AA;
5621:   PetscInt    m = *mm, n = *nn;
5622:   InsertMode  is = *isis;
5623:   Mat_SeqAIJ *a  = (Mat_SeqAIJ *)A->data;
5624:   PetscInt   *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N;
5625:   PetscInt   *imax, *ai, *ailen;
5626:   PetscInt   *aj, nonew = a->nonew, lastcol = -1;
5627:   MatScalar  *ap, value, *aa;
5628:   PetscBool   ignorezeroentries = a->ignorezeroentries;
5629:   PetscBool   roworiented       = a->roworiented;

5631:   PetscFunctionBegin;
5632:   MatCheckPreallocated(A, 1);
5633:   imax  = a->imax;
5634:   ai    = a->i;
5635:   ailen = a->ilen;
5636:   aj    = a->j;
5637:   aa    = a->a;

5639:   for (k = 0; k < m; k++) { /* loop over added rows */
5640:     row = im[k];
5641:     if (row < 0) continue;
5642:     PetscCheck(row < A->rmap->n, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Row too large");
5643:     rp   = aj + ai[row];
5644:     ap   = aa + ai[row];
5645:     rmax = imax[row];
5646:     nrow = ailen[row];
5647:     low  = 0;
5648:     high = nrow;
5649:     for (l = 0; l < n; l++) { /* loop over added columns */
5650:       if (in[l] < 0) continue;
5651:       PetscCheck(in[l] < A->cmap->n, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Column too large");
5652:       col = in[l];
5653:       if (roworiented) value = v[l + k * n];
5654:       else value = v[k + l * m];

5656:       if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;

5658:       if (col <= lastcol) low = 0;
5659:       else high = nrow;
5660:       lastcol = col;
5661:       while (high - low > 5) {
5662:         t = (low + high) / 2;
5663:         if (rp[t] > col) high = t;
5664:         else low = t;
5665:       }
5666:       for (i = low; i < high; i++) {
5667:         if (rp[i] > col) break;
5668:         if (rp[i] == col) {
5669:           if (is == ADD_VALUES) ap[i] += value;
5670:           else ap[i] = value;
5671:           goto noinsert;
5672:         }
5673:       }
5674:       if (value == 0.0 && ignorezeroentries) goto noinsert;
5675:       if (nonew == 1) goto noinsert;
5676:       PetscCheck(nonew != -1, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero in the matrix");
5677:       MatSeqXAIJReallocateAIJ(A, A->rmap->n, 1, nrow, row, col, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar);
5678:       N = nrow++ - 1;
5679:       a->nz++;
5680:       high++;
5681:       /* shift up all the later entries in this row */
5682:       for (ii = N; ii >= i; ii--) {
5683:         rp[ii + 1] = rp[ii];
5684:         ap[ii + 1] = ap[ii];
5685:       }
5686:       rp[i] = col;
5687:       ap[i] = value;
5688:     noinsert:;
5689:       low = i + 1;
5690:     }
5691:     ailen[row] = nrow;
5692:   }
5693:   PetscFunctionReturnVoid();
5694: }
5695: /* Undefining these here since they were redefined from their original definition above! No
5696:  * other PETSc functions should be defined past this point, as it is impossible to recover the
5697:  * original definitions */
5698: #undef PetscCall
5699: #undef SETERRQ