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:   if (PetscDefined(USE_DEBUG)) {
667:     PetscInt mnz = 0;

669:     for (i = 0; i < m; i++) mnz += rowlens[i];
670:     PetscCheck(nz == mnz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Row lens %" PetscInt_FMT " do not sum to nz %" PetscInt_FMT, mnz, nz);
671:   }
672:   PetscCall(PetscViewerBinaryWrite(viewer, rowlens, m, PETSC_INT));
673:   PetscCall(PetscFree(rowlens));
674:   /* store column indices */
675:   PetscCall(PetscViewerBinaryWrite(viewer, A->j, nz, PETSC_INT));
676:   /* store nonzero values */
677:   PetscCall(MatSeqAIJGetArrayRead(mat, &av));
678:   PetscCall(PetscViewerBinaryWrite(viewer, av, nz, PETSC_SCALAR));
679:   PetscCall(MatSeqAIJRestoreArrayRead(mat, &av));

681:   /* write block size option to the viewer's .info file */
682:   PetscCall(MatView_Binary_BlockSizes(mat, viewer));
683:   PetscFunctionReturn(PETSC_SUCCESS);
684: }

686: static PetscErrorCode MatView_SeqAIJ_ASCII_structonly(Mat A, PetscViewer viewer)
687: {
688:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
689:   PetscInt    i, k, m = A->rmap->N;

691:   PetscFunctionBegin;
692:   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
693:   for (i = 0; i < m; i++) {
694:     PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
695:     for (k = a->i[i]; k < a->i[i + 1]; k++) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ") ", a->j[k]));
696:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
697:   }
698:   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
699:   PetscFunctionReturn(PETSC_SUCCESS);
700: }

702: extern PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat, PetscViewer);

704: static PetscErrorCode MatView_SeqAIJ_ASCII(Mat A, PetscViewer viewer)
705: {
706:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
707:   const PetscScalar *av;
708:   PetscInt           i, j, m = A->rmap->n;
709:   const char        *name;
710:   PetscViewerFormat  format;

712:   PetscFunctionBegin;
713:   if (A->structure_only) {
714:     PetscCall(MatView_SeqAIJ_ASCII_structonly(A, viewer));
715:     PetscFunctionReturn(PETSC_SUCCESS);
716:   }

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

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

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

841:     for (i = 0; i < a->i[m]; i++) {
842:       if (PetscImaginaryPart(a->a[i]) != 0.0) {
843:         realonly = PETSC_FALSE;
844:         break;
845:       }
846:     }
847: #endif

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

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

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

977:   PetscFunctionBegin;
978:   PetscCall(PetscObjectQuery((PetscObject)A, "Zoomviewer", (PetscObject *)&viewer));
979:   PetscCall(PetscViewerGetFormat(viewer, &format));
980:   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));

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

1028:     for (i = 0; i < nz; i++) {
1029:       if (PetscAbsScalar(aa[i]) > maxv) maxv = PetscAbsScalar(aa[i]);
1030:     }
1031:     if (minv >= maxv) maxv = minv + PETSC_SMALL;
1032:     PetscCall(PetscDrawGetPopup(draw, &popup));
1033:     PetscCall(PetscDrawScalePopup(popup, minv, maxv));

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

1053: #include <petscdraw.h>
1054: static PetscErrorCode MatView_SeqAIJ_Draw(Mat A, PetscViewer viewer)
1055: {
1056:   PetscDraw draw;
1057:   PetscReal xr, yr, xl, yl, h, w;
1058:   PetscBool isnull;

1060:   PetscFunctionBegin;
1061:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1062:   PetscCall(PetscDrawIsNull(draw, &isnull));
1063:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);

1065:   xr = A->cmap->n;
1066:   yr = A->rmap->n;
1067:   h  = yr / 10.0;
1068:   w  = xr / 10.0;
1069:   xr += w;
1070:   yr += h;
1071:   xl = -w;
1072:   yl = -h;
1073:   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
1074:   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", (PetscObject)viewer));
1075:   PetscCall(PetscDrawZoom(draw, MatView_SeqAIJ_Draw_Zoom, A));
1076:   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", NULL));
1077:   PetscCall(PetscDrawSave(draw));
1078:   PetscFunctionReturn(PETSC_SUCCESS);
1079: }

1081: PetscErrorCode MatView_SeqAIJ(Mat A, PetscViewer viewer)
1082: {
1083:   PetscBool iascii, isbinary, isdraw;

1085:   PetscFunctionBegin;
1086:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1087:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1088:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1089:   if (iascii) PetscCall(MatView_SeqAIJ_ASCII(A, viewer));
1090:   else if (isbinary) PetscCall(MatView_SeqAIJ_Binary(A, viewer));
1091:   else if (isdraw) PetscCall(MatView_SeqAIJ_Draw(A, viewer));
1092:   PetscCall(MatView_SeqAIJ_Inode(A, viewer));
1093:   PetscFunctionReturn(PETSC_SUCCESS);
1094: }

1096: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A, MatAssemblyType mode)
1097: {
1098:   Mat_SeqAIJ *a      = (Mat_SeqAIJ *)A->data;
1099:   PetscInt    fshift = 0, i, *ai = a->i, *aj = a->j, *imax = a->imax;
1100:   PetscInt    m = A->rmap->n, *ip, N, *ailen = a->ilen, rmax = 0, n;
1101:   MatScalar  *aa    = a->a, *ap;
1102:   PetscReal   ratio = 0.6;

1104:   PetscFunctionBegin;
1105:   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
1106:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1107:   if (A->was_assembled && A->ass_nonzerostate == A->nonzerostate) {
1108:     /* we need to respect users asking to use or not the inodes routine in between matrix assemblies */
1109:     PetscCall(MatAssemblyEnd_SeqAIJ_Inode(A, mode));
1110:     PetscFunctionReturn(PETSC_SUCCESS);
1111:   }

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

1157:   A->info.mallocs += a->reallocs;
1158:   a->reallocs         = 0;
1159:   A->info.nz_unneeded = (PetscReal)fshift;
1160:   a->rmax             = rmax;

1162:   if (!A->structure_only) PetscCall(MatCheckCompressedRow(A, a->nonzerorowcnt, &a->compressedrow, a->i, m, ratio));
1163:   PetscCall(MatAssemblyEnd_SeqAIJ_Inode(A, mode));
1164:   PetscFunctionReturn(PETSC_SUCCESS);
1165: }

1167: static PetscErrorCode MatRealPart_SeqAIJ(Mat A)
1168: {
1169:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1170:   PetscInt    i, nz = a->nz;
1171:   MatScalar  *aa;

1173:   PetscFunctionBegin;
1174:   PetscCall(MatSeqAIJGetArray(A, &aa));
1175:   for (i = 0; i < nz; i++) aa[i] = PetscRealPart(aa[i]);
1176:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
1177:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1178:   PetscFunctionReturn(PETSC_SUCCESS);
1179: }

1181: static PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
1182: {
1183:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1184:   PetscInt    i, nz = a->nz;
1185:   MatScalar  *aa;

1187:   PetscFunctionBegin;
1188:   PetscCall(MatSeqAIJGetArray(A, &aa));
1189:   for (i = 0; i < nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1190:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
1191:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1192:   PetscFunctionReturn(PETSC_SUCCESS);
1193: }

1195: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
1196: {
1197:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1198:   MatScalar  *aa;

1200:   PetscFunctionBegin;
1201:   PetscCall(MatSeqAIJGetArrayWrite(A, &aa));
1202:   PetscCall(PetscArrayzero(aa, a->i[A->rmap->n]));
1203:   PetscCall(MatSeqAIJRestoreArrayWrite(A, &aa));
1204:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
1205:   PetscFunctionReturn(PETSC_SUCCESS);
1206: }

1208: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
1209: {
1210:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

1212:   PetscFunctionBegin;
1213:   if (A->hash_active) {
1214:     A->ops[0] = a->cops;
1215:     PetscCall(PetscHMapIJVDestroy(&a->ht));
1216:     PetscCall(PetscFree(a->dnz));
1217:     A->hash_active = PETSC_FALSE;
1218:   }

1220:   PetscCall(PetscLogObjectState((PetscObject)A, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, A->rmap->n, A->cmap->n, a->nz));
1221:   PetscCall(MatSeqXAIJFreeAIJ(A, &a->a, &a->j, &a->i));
1222:   PetscCall(ISDestroy(&a->row));
1223:   PetscCall(ISDestroy(&a->col));
1224:   PetscCall(PetscFree(a->diag));
1225:   PetscCall(PetscFree(a->ibdiag));
1226:   PetscCall(PetscFree(a->imax));
1227:   PetscCall(PetscFree(a->ilen));
1228:   PetscCall(PetscFree(a->ipre));
1229:   PetscCall(PetscFree3(a->idiag, a->mdiag, a->ssor_work));
1230:   PetscCall(PetscFree(a->solve_work));
1231:   PetscCall(ISDestroy(&a->icol));
1232:   PetscCall(PetscFree(a->saved_values));
1233:   PetscCall(PetscFree2(a->compressedrow.i, a->compressedrow.rindex));
1234:   PetscCall(MatDestroy_SeqAIJ_Inode(A));
1235:   PetscCall(PetscFree(A->data));

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

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

1306: PetscErrorCode MatSetOption_SeqAIJ(Mat A, MatOption op, PetscBool flg)
1307: {
1308:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

1310:   PetscFunctionBegin;
1311:   switch (op) {
1312:   case MAT_ROW_ORIENTED:
1313:     a->roworiented = flg;
1314:     break;
1315:   case MAT_KEEP_NONZERO_PATTERN:
1316:     a->keepnonzeropattern = flg;
1317:     break;
1318:   case MAT_NEW_NONZERO_LOCATIONS:
1319:     a->nonew = (flg ? 0 : 1);
1320:     break;
1321:   case MAT_NEW_NONZERO_LOCATION_ERR:
1322:     a->nonew = (flg ? -1 : 0);
1323:     break;
1324:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1325:     a->nonew = (flg ? -2 : 0);
1326:     break;
1327:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1328:     a->nounused = (flg ? -1 : 0);
1329:     break;
1330:   case MAT_IGNORE_ZERO_ENTRIES:
1331:     a->ignorezeroentries = flg;
1332:     break;
1333:   case MAT_SPD:
1334:   case MAT_SYMMETRIC:
1335:   case MAT_STRUCTURALLY_SYMMETRIC:
1336:   case MAT_HERMITIAN:
1337:   case MAT_SYMMETRY_ETERNAL:
1338:   case MAT_STRUCTURE_ONLY:
1339:   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
1340:   case MAT_SPD_ETERNAL:
1341:     break;
1342:   case MAT_FORCE_DIAGONAL_ENTRIES:
1343:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1344:   case MAT_USE_HASH_TABLE:
1345:     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
1346:     break;
1347:   case MAT_USE_INODES:
1348:     PetscCall(MatSetOption_SeqAIJ_Inode(A, MAT_USE_INODES, flg));
1349:     break;
1350:   case MAT_SUBMAT_SINGLEIS:
1351:     A->submat_singleis = flg;
1352:     break;
1353:   case MAT_SORTED_FULL:
1354:     if (flg) A->ops->setvalues = MatSetValues_SeqAIJ_SortedFull;
1355:     else A->ops->setvalues = MatSetValues_SeqAIJ;
1356:     break;
1357:   case MAT_FORM_EXPLICIT_TRANSPOSE:
1358:     A->form_explicit_transpose = flg;
1359:     break;
1360:   default:
1361:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
1362:   }
1363:   PetscFunctionReturn(PETSC_SUCCESS);
1364: }

1366: static PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A, Vec v)
1367: {
1368:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1369:   PetscInt           i, j, n, *ai = a->i, *aj = a->j;
1370:   PetscScalar       *x;
1371:   const PetscScalar *aa;

1373:   PetscFunctionBegin;
1374:   PetscCall(VecGetLocalSize(v, &n));
1375:   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
1376:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
1377:   if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) {
1378:     PetscInt *diag = a->diag;
1379:     PetscCall(VecGetArrayWrite(v, &x));
1380:     for (i = 0; i < n; i++) x[i] = 1.0 / aa[diag[i]];
1381:     PetscCall(VecRestoreArrayWrite(v, &x));
1382:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1383:     PetscFunctionReturn(PETSC_SUCCESS);
1384:   }

1386:   PetscCall(VecGetArrayWrite(v, &x));
1387:   for (i = 0; i < n; i++) {
1388:     x[i] = 0.0;
1389:     for (j = ai[i]; j < ai[i + 1]; j++) {
1390:       if (aj[j] == i) {
1391:         x[i] = aa[j];
1392:         break;
1393:       }
1394:     }
1395:   }
1396:   PetscCall(VecRestoreArrayWrite(v, &x));
1397:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1398:   PetscFunctionReturn(PETSC_SUCCESS);
1399: }

1401: #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>
1402: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A, Vec xx, Vec zz, Vec yy)
1403: {
1404:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1405:   const MatScalar   *aa;
1406:   PetscScalar       *y;
1407:   const PetscScalar *x;
1408:   PetscInt           m = A->rmap->n;
1409: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1410:   const MatScalar  *v;
1411:   PetscScalar       alpha;
1412:   PetscInt          n, i, j;
1413:   const PetscInt   *idx, *ii, *ridx = NULL;
1414:   Mat_CompressedRow cprow    = a->compressedrow;
1415:   PetscBool         usecprow = cprow.use;
1416: #endif

1418:   PetscFunctionBegin;
1419:   if (zz != yy) PetscCall(VecCopy(zz, yy));
1420:   PetscCall(VecGetArrayRead(xx, &x));
1421:   PetscCall(VecGetArray(yy, &y));
1422:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));

1424: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1425:   fortranmulttransposeaddaij_(&m, x, a->i, a->j, aa, y);
1426: #else
1427:   if (usecprow) {
1428:     m    = cprow.nrows;
1429:     ii   = cprow.i;
1430:     ridx = cprow.rindex;
1431:   } else {
1432:     ii = a->i;
1433:   }
1434:   for (i = 0; i < m; i++) {
1435:     idx = a->j + ii[i];
1436:     v   = aa + ii[i];
1437:     n   = ii[i + 1] - ii[i];
1438:     if (usecprow) {
1439:       alpha = x[ridx[i]];
1440:     } else {
1441:       alpha = x[i];
1442:     }
1443:     for (j = 0; j < n; j++) y[idx[j]] += alpha * v[j];
1444:   }
1445: #endif
1446:   PetscCall(PetscLogFlops(2.0 * a->nz));
1447:   PetscCall(VecRestoreArrayRead(xx, &x));
1448:   PetscCall(VecRestoreArray(yy, &y));
1449:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1450:   PetscFunctionReturn(PETSC_SUCCESS);
1451: }

1453: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A, Vec xx, Vec yy)
1454: {
1455:   PetscFunctionBegin;
1456:   PetscCall(VecSet(yy, 0.0));
1457:   PetscCall(MatMultTransposeAdd_SeqAIJ(A, xx, yy, yy));
1458:   PetscFunctionReturn(PETSC_SUCCESS);
1459: }

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

1463: PetscErrorCode MatMult_SeqAIJ(Mat A, Vec xx, Vec yy)
1464: {
1465:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1466:   PetscScalar       *y;
1467:   const PetscScalar *x;
1468:   const MatScalar   *a_a;
1469:   PetscInt           m = A->rmap->n;
1470:   const PetscInt    *ii, *ridx = NULL;
1471:   PetscBool          usecprow = a->compressedrow.use;

1473: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1474:   #pragma disjoint(*x, *y, *aa)
1475: #endif

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

1523: // HACK!!!!! Used by src/mat/tests/ex170.c
1524: PETSC_EXTERN PetscErrorCode MatMultMax_SeqAIJ(Mat A, Vec xx, Vec yy)
1525: {
1526:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1527:   PetscScalar       *y;
1528:   const PetscScalar *x;
1529:   const MatScalar   *aa, *a_a;
1530:   PetscInt           m = A->rmap->n;
1531:   const PetscInt    *aj, *ii, *ridx   = NULL;
1532:   PetscInt           n, i, nonzerorow = 0;
1533:   PetscScalar        sum;
1534:   PetscBool          usecprow = a->compressedrow.use;

1536: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1537:   #pragma disjoint(*x, *y, *aa)
1538: #endif

1540:   PetscFunctionBegin;
1541:   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1542:   PetscCall(VecGetArrayRead(xx, &x));
1543:   PetscCall(VecGetArray(yy, &y));
1544:   if (usecprow) { /* use compressed row format */
1545:     m    = a->compressedrow.nrows;
1546:     ii   = a->compressedrow.i;
1547:     ridx = a->compressedrow.rindex;
1548:     for (i = 0; i < m; i++) {
1549:       n   = ii[i + 1] - ii[i];
1550:       aj  = a->j + ii[i];
1551:       aa  = a_a + ii[i];
1552:       sum = 0.0;
1553:       nonzerorow += (n > 0);
1554:       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1555:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1556:       y[*ridx++] = sum;
1557:     }
1558:   } else { /* do not use compressed row format */
1559:     ii = a->i;
1560:     for (i = 0; i < m; i++) {
1561:       n   = ii[i + 1] - ii[i];
1562:       aj  = a->j + ii[i];
1563:       aa  = a_a + ii[i];
1564:       sum = 0.0;
1565:       nonzerorow += (n > 0);
1566:       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1567:       y[i] = sum;
1568:     }
1569:   }
1570:   PetscCall(PetscLogFlops(2.0 * a->nz - nonzerorow));
1571:   PetscCall(VecRestoreArrayRead(xx, &x));
1572:   PetscCall(VecRestoreArray(yy, &y));
1573:   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1574:   PetscFunctionReturn(PETSC_SUCCESS);
1575: }

1577: // HACK!!!!! Used by src/mat/tests/ex170.c
1578: PETSC_EXTERN PetscErrorCode MatMultAddMax_SeqAIJ(Mat A, Vec xx, Vec yy, Vec zz)
1579: {
1580:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1581:   PetscScalar       *y, *z;
1582:   const PetscScalar *x;
1583:   const MatScalar   *aa, *a_a;
1584:   PetscInt           m = A->rmap->n, *aj, *ii;
1585:   PetscInt           n, i, *ridx = NULL;
1586:   PetscScalar        sum;
1587:   PetscBool          usecprow = a->compressedrow.use;

1589:   PetscFunctionBegin;
1590:   PetscCall(MatSeqAIJGetArrayRead(A, &a_a));
1591:   PetscCall(VecGetArrayRead(xx, &x));
1592:   PetscCall(VecGetArrayPair(yy, zz, &y, &z));
1593:   if (usecprow) { /* use compressed row format */
1594:     if (zz != yy) PetscCall(PetscArraycpy(z, y, m));
1595:     m    = a->compressedrow.nrows;
1596:     ii   = a->compressedrow.i;
1597:     ridx = a->compressedrow.rindex;
1598:     for (i = 0; i < m; i++) {
1599:       n   = ii[i + 1] - ii[i];
1600:       aj  = a->j + ii[i];
1601:       aa  = a_a + ii[i];
1602:       sum = y[*ridx];
1603:       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1604:       z[*ridx++] = sum;
1605:     }
1606:   } else { /* do not use compressed row format */
1607:     ii = a->i;
1608:     for (i = 0; i < m; i++) {
1609:       n   = ii[i + 1] - ii[i];
1610:       aj  = a->j + ii[i];
1611:       aa  = a_a + ii[i];
1612:       sum = y[i];
1613:       PetscSparseDenseMaxDot(sum, x, aa, aj, n);
1614:       z[i] = sum;
1615:     }
1616:   }
1617:   PetscCall(PetscLogFlops(2.0 * a->nz));
1618:   PetscCall(VecRestoreArrayRead(xx, &x));
1619:   PetscCall(VecRestoreArrayPair(yy, zz, &y, &z));
1620:   PetscCall(MatSeqAIJRestoreArrayRead(A, &a_a));
1621:   PetscFunctionReturn(PETSC_SUCCESS);
1622: }

1624: #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h>
1625: PetscErrorCode MatMultAdd_SeqAIJ(Mat A, Vec xx, Vec yy, Vec zz)
1626: {
1627:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1628:   PetscScalar       *y, *z;
1629:   const PetscScalar *x;
1630:   const MatScalar   *a_a;
1631:   const PetscInt    *ii, *ridx = NULL;
1632:   PetscInt           m        = A->rmap->n;
1633:   PetscBool          usecprow = a->compressedrow.use;

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

1679: /*
1680:      Adds diagonal pointers to sparse matrix nonzero structure.
1681: */
1682: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1683: {
1684:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1685:   PetscInt    i, j, m = A->rmap->n;
1686:   PetscBool   alreadySet = PETSC_TRUE;

1688:   PetscFunctionBegin;
1689:   if (!a->diag) {
1690:     PetscCall(PetscMalloc1(m, &a->diag));
1691:     alreadySet = PETSC_FALSE;
1692:   }
1693:   for (i = 0; i < A->rmap->n; i++) {
1694:     /* If A's diagonal is already correctly set, this fast track enables cheap and repeated MatMarkDiagonal_SeqAIJ() calls */
1695:     if (alreadySet) {
1696:       PetscInt pos = a->diag[i];
1697:       if (pos >= a->i[i] && pos < a->i[i + 1] && a->j[pos] == i) continue;
1698:     }

1700:     a->diag[i] = a->i[i + 1];
1701:     for (j = a->i[i]; j < a->i[i + 1]; j++) {
1702:       if (a->j[j] == i) {
1703:         a->diag[i] = j;
1704:         break;
1705:       }
1706:     }
1707:   }
1708:   PetscFunctionReturn(PETSC_SUCCESS);
1709: }

1711: static PetscErrorCode MatShift_SeqAIJ(Mat A, PetscScalar v)
1712: {
1713:   Mat_SeqAIJ     *a    = (Mat_SeqAIJ *)A->data;
1714:   const PetscInt *diag = (const PetscInt *)a->diag;
1715:   const PetscInt *ii   = (const PetscInt *)a->i;
1716:   PetscInt        i, *mdiag = NULL;
1717:   PetscInt        cnt = 0; /* how many diagonals are missing */

1719:   PetscFunctionBegin;
1720:   if (!A->preallocated || !a->nz) {
1721:     PetscCall(MatSeqAIJSetPreallocation(A, 1, NULL));
1722:     PetscCall(MatShift_Basic(A, v));
1723:     PetscFunctionReturn(PETSC_SUCCESS);
1724:   }

1726:   if (a->diagonaldense) {
1727:     cnt = 0;
1728:   } else {
1729:     PetscCall(PetscCalloc1(A->rmap->n, &mdiag));
1730:     for (i = 0; i < A->rmap->n; i++) {
1731:       if (i < A->cmap->n && diag[i] >= ii[i + 1]) { /* 'out of range' rows never have diagonals */
1732:         cnt++;
1733:         mdiag[i] = 1;
1734:       }
1735:     }
1736:   }
1737:   if (!cnt) {
1738:     PetscCall(MatShift_Basic(A, v));
1739:   } else {
1740:     PetscScalar       *olda = a->a; /* preserve pointers to current matrix nonzeros structure and values */
1741:     PetscInt          *oldj = a->j, *oldi = a->i;
1742:     PetscBool          free_a = a->free_a, free_ij = a->free_ij;
1743:     const PetscScalar *Aa;

1745:     PetscCall(MatSeqAIJGetArrayRead(A, &Aa)); // sync the host
1746:     PetscCall(MatSeqAIJRestoreArrayRead(A, &Aa));

1748:     a->a = NULL;
1749:     a->j = NULL;
1750:     a->i = NULL;
1751:     /* increase the values in imax for each row where a diagonal is being inserted then reallocate the matrix data structures */
1752:     for (i = 0; i < PetscMin(A->rmap->n, A->cmap->n); i++) a->imax[i] += mdiag[i];
1753:     PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(A, 0, a->imax));

1755:     /* copy old values into new matrix data structure */
1756:     for (i = 0; i < A->rmap->n; i++) {
1757:       PetscCall(MatSetValues(A, 1, &i, a->imax[i] - mdiag[i], &oldj[oldi[i]], &olda[oldi[i]], ADD_VALUES));
1758:       if (i < A->cmap->n) PetscCall(MatSetValue(A, i, i, v, ADD_VALUES));
1759:     }
1760:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1761:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1762:     if (free_a) PetscCall(PetscShmgetDeallocateArray((void **)&olda));
1763:     if (free_ij) PetscCall(PetscShmgetDeallocateArray((void **)&oldj));
1764:     if (free_ij) PetscCall(PetscShmgetDeallocateArray((void **)&oldi));
1765:   }
1766:   PetscCall(PetscFree(mdiag));
1767:   a->diagonaldense = PETSC_TRUE;
1768:   PetscFunctionReturn(PETSC_SUCCESS);
1769: }

1771: /*
1772:      Checks for missing diagonals
1773: */
1774: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A, PetscBool *missing, PetscInt *d)
1775: {
1776:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
1777:   PetscInt   *diag, *ii = a->i, i;

1779:   PetscFunctionBegin;
1780:   *missing = PETSC_FALSE;
1781:   if (A->rmap->n > 0 && !ii) {
1782:     *missing = PETSC_TRUE;
1783:     if (d) *d = 0;
1784:     PetscCall(PetscInfo(A, "Matrix has no entries therefore is missing diagonal\n"));
1785:   } else {
1786:     PetscInt n;
1787:     n    = PetscMin(A->rmap->n, A->cmap->n);
1788:     diag = a->diag;
1789:     for (i = 0; i < n; i++) {
1790:       if (diag[i] >= ii[i + 1]) {
1791:         *missing = PETSC_TRUE;
1792:         if (d) *d = i;
1793:         PetscCall(PetscInfo(A, "Matrix is missing diagonal number %" PetscInt_FMT "\n", i));
1794:         break;
1795:       }
1796:     }
1797:   }
1798:   PetscFunctionReturn(PETSC_SUCCESS);
1799: }

1801: #include <petscblaslapack.h>
1802: #include <petsc/private/kernels/blockinvert.h>

1804: /*
1805:     Note that values is allocated externally by the PC and then passed into this routine
1806: */
1807: static PetscErrorCode MatInvertVariableBlockDiagonal_SeqAIJ(Mat A, PetscInt nblocks, const PetscInt *bsizes, PetscScalar *diag)
1808: {
1809:   PetscInt        n = A->rmap->n, i, ncnt = 0, *indx, j, bsizemax = 0, *v_pivots;
1810:   PetscBool       allowzeropivot, zeropivotdetected = PETSC_FALSE;
1811:   const PetscReal shift = 0.0;
1812:   PetscInt        ipvt[5];
1813:   PetscCount      flops = 0;
1814:   PetscScalar     work[25], *v_work;

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

1876: /*
1877:    Negative shift indicates do not generate an error if there is a zero diagonal, just invert it anyways
1878: */
1879: static PetscErrorCode MatInvertDiagonal_SeqAIJ(Mat A, PetscScalar omega, PetscScalar fshift)
1880: {
1881:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
1882:   PetscInt         i, *diag, m = A->rmap->n;
1883:   const MatScalar *v;
1884:   PetscScalar     *idiag, *mdiag;

1886:   PetscFunctionBegin;
1887:   if (a->idiagvalid) PetscFunctionReturn(PETSC_SUCCESS);
1888:   PetscCall(MatMarkDiagonal_SeqAIJ(A));
1889:   diag = a->diag;
1890:   if (!a->idiag) { PetscCall(PetscMalloc3(m, &a->idiag, m, &a->mdiag, m, &a->ssor_work)); }

1892:   mdiag = a->mdiag;
1893:   idiag = a->idiag;
1894:   PetscCall(MatSeqAIJGetArrayRead(A, &v));
1895:   if (omega == 1.0 && PetscRealPart(fshift) <= 0.0) {
1896:     for (i = 0; i < m; i++) {
1897:       mdiag[i] = v[diag[i]];
1898:       if (!PetscAbsScalar(mdiag[i])) { /* zero diagonal */
1899:         if (PetscRealPart(fshift)) {
1900:           PetscCall(PetscInfo(A, "Zero diagonal on row %" PetscInt_FMT "\n", i));
1901:           A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1902:           A->factorerror_zeropivot_value = 0.0;
1903:           A->factorerror_zeropivot_row   = i;
1904:         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Zero diagonal on row %" PetscInt_FMT, i);
1905:       }
1906:       idiag[i] = 1.0 / v[diag[i]];
1907:     }
1908:     PetscCall(PetscLogFlops(m));
1909:   } else {
1910:     for (i = 0; i < m; i++) {
1911:       mdiag[i] = v[diag[i]];
1912:       idiag[i] = omega / (fshift + v[diag[i]]);
1913:     }
1914:     PetscCall(PetscLogFlops(2.0 * m));
1915:   }
1916:   a->idiagvalid = PETSC_TRUE;
1917:   PetscCall(MatSeqAIJRestoreArrayRead(A, &v));
1918:   PetscFunctionReturn(PETSC_SUCCESS);
1919: }

1921: PetscErrorCode MatSOR_SeqAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
1922: {
1923:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
1924:   PetscScalar       *x, d, sum, *t, scale;
1925:   const MatScalar   *v, *idiag = NULL, *mdiag, *aa;
1926:   const PetscScalar *b, *bs, *xb, *ts;
1927:   PetscInt           n, m = A->rmap->n, i;
1928:   const PetscInt    *idx, *diag;

1930:   PetscFunctionBegin;
1931:   if (a->inode.use && a->inode.checked && omega == 1.0 && fshift == 0.0) {
1932:     PetscCall(MatSOR_SeqAIJ_Inode(A, bb, omega, flag, fshift, its, lits, xx));
1933:     PetscFunctionReturn(PETSC_SUCCESS);
1934:   }
1935:   its = its * lits;

1937:   if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
1938:   if (!a->idiagvalid) PetscCall(MatInvertDiagonal_SeqAIJ(A, omega, fshift));
1939:   a->fshift = fshift;
1940:   a->omega  = omega;

1942:   diag  = a->diag;
1943:   t     = a->ssor_work;
1944:   idiag = a->idiag;
1945:   mdiag = a->mdiag;

1947:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
1948:   PetscCall(VecGetArray(xx, &x));
1949:   PetscCall(VecGetArrayRead(bb, &b));
1950:   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1951:   if (flag == SOR_APPLY_UPPER) {
1952:     /* apply (U + D/omega) to the vector */
1953:     bs = b;
1954:     for (i = 0; i < m; i++) {
1955:       d   = fshift + mdiag[i];
1956:       n   = a->i[i + 1] - diag[i] - 1;
1957:       idx = a->j + diag[i] + 1;
1958:       v   = aa + diag[i] + 1;
1959:       sum = b[i] * d / omega;
1960:       PetscSparseDensePlusDot(sum, bs, v, idx, n);
1961:       x[i] = sum;
1962:     }
1963:     PetscCall(VecRestoreArray(xx, &x));
1964:     PetscCall(VecRestoreArrayRead(bb, &b));
1965:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
1966:     PetscCall(PetscLogFlops(a->nz));
1967:     PetscFunctionReturn(PETSC_SUCCESS);
1968:   }

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

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

1977:     to a vector efficiently using Eisenstat's trick.
1978:     */
1979:     scale = (2.0 / omega) - 1.0;

1981:     /*  x = (E + U)^{-1} b */
1982:     for (i = m - 1; i >= 0; i--) {
1983:       n   = a->i[i + 1] - diag[i] - 1;
1984:       idx = a->j + diag[i] + 1;
1985:       v   = aa + diag[i] + 1;
1986:       sum = b[i];
1987:       PetscSparseDenseMinusDot(sum, x, v, idx, n);
1988:       x[i] = sum * idiag[i];
1989:     }

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

1995:     /*  t = (E + L)^{-1}t */
1996:     ts   = t;
1997:     diag = a->diag;
1998:     for (i = 0; i < m; i++) {
1999:       n   = diag[i] - a->i[i];
2000:       idx = a->j + a->i[i];
2001:       v   = aa + a->i[i];
2002:       sum = t[i];
2003:       PetscSparseDenseMinusDot(sum, ts, v, idx, n);
2004:       t[i] = sum * idiag[i];
2005:       /*  x = x + t */
2006:       x[i] += t[i];
2007:     }

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

2096: static PetscErrorCode MatGetInfo_SeqAIJ(Mat A, MatInfoType flag, MatInfo *info)
2097: {
2098:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

2100:   PetscFunctionBegin;
2101:   info->block_size   = 1.0;
2102:   info->nz_allocated = a->maxnz;
2103:   info->nz_used      = a->nz;
2104:   info->nz_unneeded  = (a->maxnz - a->nz);
2105:   info->assemblies   = A->num_ass;
2106:   info->mallocs      = A->info.mallocs;
2107:   info->memory       = 0; /* REVIEW ME */
2108:   if (A->factortype) {
2109:     info->fill_ratio_given  = A->info.fill_ratio_given;
2110:     info->fill_ratio_needed = A->info.fill_ratio_needed;
2111:     info->factor_mallocs    = A->info.factor_mallocs;
2112:   } else {
2113:     info->fill_ratio_given  = 0;
2114:     info->fill_ratio_needed = 0;
2115:     info->factor_mallocs    = 0;
2116:   }
2117:   PetscFunctionReturn(PETSC_SUCCESS);
2118: }

2120: static PetscErrorCode MatZeroRows_SeqAIJ(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
2121: {
2122:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2123:   PetscInt           i, m = A->rmap->n - 1;
2124:   const PetscScalar *xx;
2125:   PetscScalar       *bb, *aa;
2126:   PetscInt           d = 0;

2128:   PetscFunctionBegin;
2129:   if (x && b) {
2130:     PetscCall(VecGetArrayRead(x, &xx));
2131:     PetscCall(VecGetArray(b, &bb));
2132:     for (i = 0; i < N; i++) {
2133:       PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2134:       if (rows[i] >= A->cmap->n) continue;
2135:       bb[rows[i]] = diag * xx[rows[i]];
2136:     }
2137:     PetscCall(VecRestoreArrayRead(x, &xx));
2138:     PetscCall(VecRestoreArray(b, &bb));
2139:   }

2141:   PetscCall(MatSeqAIJGetArray(A, &aa));
2142:   if (a->keepnonzeropattern) {
2143:     for (i = 0; i < N; i++) {
2144:       PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2145:       PetscCall(PetscArrayzero(&aa[a->i[rows[i]]], a->ilen[rows[i]]));
2146:     }
2147:     if (diag != 0.0) {
2148:       for (i = 0; i < N; i++) {
2149:         d = rows[i];
2150:         if (rows[i] >= A->cmap->n) continue;
2151:         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);
2152:       }
2153:       for (i = 0; i < N; i++) {
2154:         if (rows[i] >= A->cmap->n) continue;
2155:         aa[a->diag[rows[i]]] = diag;
2156:       }
2157:     }
2158:   } else {
2159:     if (diag != 0.0) {
2160:       for (i = 0; i < N; i++) {
2161:         PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2162:         if (a->ilen[rows[i]] > 0) {
2163:           if (rows[i] >= A->cmap->n) {
2164:             a->ilen[rows[i]] = 0;
2165:           } else {
2166:             a->ilen[rows[i]]    = 1;
2167:             aa[a->i[rows[i]]]   = diag;
2168:             a->j[a->i[rows[i]]] = rows[i];
2169:           }
2170:         } else if (rows[i] < A->cmap->n) { /* in case row was completely empty */
2171:           PetscCall(MatSetValues_SeqAIJ(A, 1, &rows[i], 1, &rows[i], &diag, INSERT_VALUES));
2172:         }
2173:       }
2174:     } else {
2175:       for (i = 0; i < N; i++) {
2176:         PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2177:         a->ilen[rows[i]] = 0;
2178:       }
2179:     }
2180:     A->nonzerostate++;
2181:   }
2182:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
2183:   PetscUseTypeMethod(A, assemblyend, MAT_FINAL_ASSEMBLY);
2184:   PetscFunctionReturn(PETSC_SUCCESS);
2185: }

2187: static PetscErrorCode MatZeroRowsColumns_SeqAIJ(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
2188: {
2189:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2190:   PetscInt           i, j, m = A->rmap->n - 1, d = 0;
2191:   PetscBool          missing, *zeroed, vecs = PETSC_FALSE;
2192:   const PetscScalar *xx;
2193:   PetscScalar       *bb, *aa;

2195:   PetscFunctionBegin;
2196:   if (!N) PetscFunctionReturn(PETSC_SUCCESS);
2197:   PetscCall(MatSeqAIJGetArray(A, &aa));
2198:   if (x && b) {
2199:     PetscCall(VecGetArrayRead(x, &xx));
2200:     PetscCall(VecGetArray(b, &bb));
2201:     vecs = PETSC_TRUE;
2202:   }
2203:   PetscCall(PetscCalloc1(A->rmap->n, &zeroed));
2204:   for (i = 0; i < N; i++) {
2205:     PetscCheck(rows[i] >= 0 && rows[i] <= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", rows[i]);
2206:     PetscCall(PetscArrayzero(PetscSafePointerPlusOffset(aa, a->i[rows[i]]), a->ilen[rows[i]]));

2208:     zeroed[rows[i]] = PETSC_TRUE;
2209:   }
2210:   for (i = 0; i < A->rmap->n; i++) {
2211:     if (!zeroed[i]) {
2212:       for (j = a->i[i]; j < a->i[i + 1]; j++) {
2213:         if (a->j[j] < A->rmap->n && zeroed[a->j[j]]) {
2214:           if (vecs) bb[i] -= aa[j] * xx[a->j[j]];
2215:           aa[j] = 0.0;
2216:         }
2217:       }
2218:     } else if (vecs && i < A->cmap->N) bb[i] = diag * xx[i];
2219:   }
2220:   if (x && b) {
2221:     PetscCall(VecRestoreArrayRead(x, &xx));
2222:     PetscCall(VecRestoreArray(b, &bb));
2223:   }
2224:   PetscCall(PetscFree(zeroed));
2225:   if (diag != 0.0) {
2226:     PetscCall(MatMissingDiagonal_SeqAIJ(A, &missing, &d));
2227:     if (missing) {
2228:       for (i = 0; i < N; i++) {
2229:         if (rows[i] >= A->cmap->N) continue;
2230:         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]);
2231:         PetscCall(MatSetValues_SeqAIJ(A, 1, &rows[i], 1, &rows[i], &diag, INSERT_VALUES));
2232:       }
2233:     } else {
2234:       for (i = 0; i < N; i++) aa[a->diag[rows[i]]] = diag;
2235:     }
2236:   }
2237:   PetscCall(MatSeqAIJRestoreArray(A, &aa));
2238:   PetscUseTypeMethod(A, assemblyend, MAT_FINAL_ASSEMBLY);
2239:   PetscFunctionReturn(PETSC_SUCCESS);
2240: }

2242: PetscErrorCode MatGetRow_SeqAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
2243: {
2244:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2245:   const PetscScalar *aa;

2247:   PetscFunctionBegin;
2248:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2249:   *nz = a->i[row + 1] - a->i[row];
2250:   if (v) *v = PetscSafePointerPlusOffset((PetscScalar *)aa, a->i[row]);
2251:   if (idx) {
2252:     if (*nz && a->j) *idx = a->j + a->i[row];
2253:     else *idx = NULL;
2254:   }
2255:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2256:   PetscFunctionReturn(PETSC_SUCCESS);
2257: }

2259: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
2260: {
2261:   PetscFunctionBegin;
2262:   PetscFunctionReturn(PETSC_SUCCESS);
2263: }

2265: static PetscErrorCode MatNorm_SeqAIJ(Mat A, NormType type, PetscReal *nrm)
2266: {
2267:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
2268:   const MatScalar *v;
2269:   PetscReal        sum = 0.0;
2270:   PetscInt         i, j;

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

2317: static PetscErrorCode MatIsTranspose_SeqAIJ(Mat A, Mat B, PetscReal tol, PetscBool *f)
2318: {
2319:   Mat_SeqAIJ      *aij = (Mat_SeqAIJ *)A->data, *bij = (Mat_SeqAIJ *)B->data;
2320:   PetscInt        *adx, *bdx, *aii, *bii, *aptr, *bptr;
2321:   const MatScalar *va, *vb;
2322:   PetscInt         ma, na, mb, nb, i;

2324:   PetscFunctionBegin;
2325:   PetscCall(MatGetSize(A, &ma, &na));
2326:   PetscCall(MatGetSize(B, &mb, &nb));
2327:   if (ma != nb || na != mb) {
2328:     *f = PETSC_FALSE;
2329:     PetscFunctionReturn(PETSC_SUCCESS);
2330:   }
2331:   PetscCall(MatSeqAIJGetArrayRead(A, &va));
2332:   PetscCall(MatSeqAIJGetArrayRead(B, &vb));
2333:   aii = aij->i;
2334:   bii = bij->i;
2335:   adx = aij->j;
2336:   bdx = bij->j;
2337:   PetscCall(PetscMalloc1(ma, &aptr));
2338:   PetscCall(PetscMalloc1(mb, &bptr));
2339:   for (i = 0; i < ma; i++) aptr[i] = aii[i];
2340:   for (i = 0; i < mb; i++) bptr[i] = bii[i];

2342:   *f = PETSC_TRUE;
2343:   for (i = 0; i < ma; i++) {
2344:     while (aptr[i] < aii[i + 1]) {
2345:       PetscInt    idc, idr;
2346:       PetscScalar vc, vr;
2347:       /* column/row index/value */
2348:       idc = adx[aptr[i]];
2349:       idr = bdx[bptr[idc]];
2350:       vc  = va[aptr[i]];
2351:       vr  = vb[bptr[idc]];
2352:       if (i != idr || PetscAbsScalar(vc - vr) > tol) {
2353:         *f = PETSC_FALSE;
2354:         goto done;
2355:       } else {
2356:         aptr[i]++;
2357:         if (B || i != idc) bptr[idc]++;
2358:       }
2359:     }
2360:   }
2361: done:
2362:   PetscCall(PetscFree(aptr));
2363:   PetscCall(PetscFree(bptr));
2364:   PetscCall(MatSeqAIJRestoreArrayRead(A, &va));
2365:   PetscCall(MatSeqAIJRestoreArrayRead(B, &vb));
2366:   PetscFunctionReturn(PETSC_SUCCESS);
2367: }

2369: static PetscErrorCode MatIsHermitianTranspose_SeqAIJ(Mat A, Mat B, PetscReal tol, PetscBool *f)
2370: {
2371:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data, *bij = (Mat_SeqAIJ *)B->data;
2372:   PetscInt   *adx, *bdx, *aii, *bii, *aptr, *bptr;
2373:   MatScalar  *va, *vb;
2374:   PetscInt    ma, na, mb, nb, i;

2376:   PetscFunctionBegin;
2377:   PetscCall(MatGetSize(A, &ma, &na));
2378:   PetscCall(MatGetSize(B, &mb, &nb));
2379:   if (ma != nb || na != mb) {
2380:     *f = PETSC_FALSE;
2381:     PetscFunctionReturn(PETSC_SUCCESS);
2382:   }
2383:   aii = aij->i;
2384:   bii = bij->i;
2385:   adx = aij->j;
2386:   bdx = bij->j;
2387:   va  = aij->a;
2388:   vb  = bij->a;
2389:   PetscCall(PetscMalloc1(ma, &aptr));
2390:   PetscCall(PetscMalloc1(mb, &bptr));
2391:   for (i = 0; i < ma; i++) aptr[i] = aii[i];
2392:   for (i = 0; i < mb; i++) bptr[i] = bii[i];

2394:   *f = PETSC_TRUE;
2395:   for (i = 0; i < ma; i++) {
2396:     while (aptr[i] < aii[i + 1]) {
2397:       PetscInt    idc, idr;
2398:       PetscScalar vc, vr;
2399:       /* column/row index/value */
2400:       idc = adx[aptr[i]];
2401:       idr = bdx[bptr[idc]];
2402:       vc  = va[aptr[i]];
2403:       vr  = vb[bptr[idc]];
2404:       if (i != idr || PetscAbsScalar(vc - PetscConj(vr)) > tol) {
2405:         *f = PETSC_FALSE;
2406:         goto done;
2407:       } else {
2408:         aptr[i]++;
2409:         if (B || i != idc) bptr[idc]++;
2410:       }
2411:     }
2412:   }
2413: done:
2414:   PetscCall(PetscFree(aptr));
2415:   PetscCall(PetscFree(bptr));
2416:   PetscFunctionReturn(PETSC_SUCCESS);
2417: }

2419: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A, Vec ll, Vec rr)
2420: {
2421:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2422:   const PetscScalar *l, *r;
2423:   PetscScalar        x;
2424:   MatScalar         *v;
2425:   PetscInt           i, j, m = A->rmap->n, n = A->cmap->n, M, nz = a->nz;
2426:   const PetscInt    *jj;

2428:   PetscFunctionBegin;
2429:   if (ll) {
2430:     /* The local size is used so that VecMPI can be passed to this routine
2431:        by MatDiagonalScale_MPIAIJ */
2432:     PetscCall(VecGetLocalSize(ll, &m));
2433:     PetscCheck(m == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Left scaling vector wrong length");
2434:     PetscCall(VecGetArrayRead(ll, &l));
2435:     PetscCall(MatSeqAIJGetArray(A, &v));
2436:     for (i = 0; i < m; i++) {
2437:       x = l[i];
2438:       M = a->i[i + 1] - a->i[i];
2439:       for (j = 0; j < M; j++) (*v++) *= x;
2440:     }
2441:     PetscCall(VecRestoreArrayRead(ll, &l));
2442:     PetscCall(PetscLogFlops(nz));
2443:     PetscCall(MatSeqAIJRestoreArray(A, &v));
2444:   }
2445:   if (rr) {
2446:     PetscCall(VecGetLocalSize(rr, &n));
2447:     PetscCheck(n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Right scaling vector wrong length");
2448:     PetscCall(VecGetArrayRead(rr, &r));
2449:     PetscCall(MatSeqAIJGetArray(A, &v));
2450:     jj = a->j;
2451:     for (i = 0; i < nz; i++) (*v++) *= r[*jj++];
2452:     PetscCall(MatSeqAIJRestoreArray(A, &v));
2453:     PetscCall(VecRestoreArrayRead(rr, &r));
2454:     PetscCall(PetscLogFlops(nz));
2455:   }
2456:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
2457:   PetscFunctionReturn(PETSC_SUCCESS);
2458: }

2460: PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat A, IS isrow, IS iscol, PetscInt csize, MatReuse scall, Mat *B)
2461: {
2462:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data, *c;
2463:   PetscInt          *smap, i, k, kstart, kend, oldcols = A->cmap->n, *lens;
2464:   PetscInt           row, mat_i, *mat_j, tcol, first, step, *mat_ilen, sum, lensi;
2465:   const PetscInt    *irow, *icol;
2466:   const PetscScalar *aa;
2467:   PetscInt           nrows, ncols;
2468:   PetscInt          *starts, *j_new, *i_new, *aj = a->j, *ai = a->i, ii, *ailen = a->ilen;
2469:   MatScalar         *a_new, *mat_a, *c_a;
2470:   Mat                C;
2471:   PetscBool          stride;

2473:   PetscFunctionBegin;
2474:   PetscCall(ISGetIndices(isrow, &irow));
2475:   PetscCall(ISGetLocalSize(isrow, &nrows));
2476:   PetscCall(ISGetLocalSize(iscol, &ncols));

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

2525:     /* loop over rows inserting into submatrix */
2526:     PetscCall(MatSeqAIJGetArrayWrite(C, &a_new)); // Not 'a_new = c->a-new', since that raw usage ignores offload state of C
2527:     j_new = c->j;
2528:     i_new = c->i;
2529:     PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2530:     for (i = 0; i < nrows; i++) {
2531:       ii    = starts[i];
2532:       lensi = lens[i];
2533:       if (lensi) {
2534:         for (k = 0; k < lensi; k++) *j_new++ = aj[ii + k] - first;
2535:         PetscCall(PetscArraycpy(a_new, aa + starts[i], lensi));
2536:         a_new += lensi;
2537:       }
2538:       i_new[i + 1] = i_new[i] + lensi;
2539:       c->ilen[i]   = lensi;
2540:     }
2541:     PetscCall(MatSeqAIJRestoreArrayWrite(C, &a_new)); // Set C's offload state properly
2542:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2543:     PetscCall(PetscFree2(lens, starts));
2544:   } else {
2545:     PetscCall(ISGetIndices(iscol, &icol));
2546:     PetscCall(PetscCalloc1(oldcols, &smap));
2547:     PetscCall(PetscMalloc1(1 + nrows, &lens));
2548:     for (i = 0; i < ncols; i++) {
2549:       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);
2550:       smap[icol[i]] = i + 1;
2551:     }

2553:     /* determine lens of each row */
2554:     for (i = 0; i < nrows; i++) {
2555:       kstart  = ai[irow[i]];
2556:       kend    = kstart + a->ilen[irow[i]];
2557:       lens[i] = 0;
2558:       for (k = kstart; k < kend; k++) {
2559:         if (smap[aj[k]]) lens[i]++;
2560:       }
2561:     }
2562:     /* Create and fill new matrix */
2563:     if (scall == MAT_REUSE_MATRIX) {
2564:       PetscBool equal;

2566:       c = (Mat_SeqAIJ *)((*B)->data);
2567:       PetscCheck((*B)->rmap->n == nrows && (*B)->cmap->n == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong size");
2568:       PetscCall(PetscArraycmp(c->ilen, lens, (*B)->rmap->n, &equal));
2569:       PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong number of nonzeros");
2570:       PetscCall(PetscArrayzero(c->ilen, (*B)->rmap->n));
2571:       C = *B;
2572:     } else {
2573:       PetscInt rbs, cbs;
2574:       PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
2575:       PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE));
2576:       PetscCall(ISGetBlockSize(isrow, &rbs));
2577:       PetscCall(ISGetBlockSize(iscol, &cbs));
2578:       if (rbs > 1 || cbs > 1) PetscCall(MatSetBlockSizes(C, rbs, cbs));
2579:       PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
2580:       PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(C, 0, lens));
2581:     }
2582:     PetscCall(MatSeqAIJGetArrayRead(A, &aa));

2584:     c = (Mat_SeqAIJ *)C->data;
2585:     PetscCall(MatSeqAIJGetArrayWrite(C, &c_a)); // Not 'c->a', since that raw usage ignores offload state of C
2586:     for (i = 0; i < nrows; i++) {
2587:       row      = irow[i];
2588:       kstart   = ai[row];
2589:       kend     = kstart + a->ilen[row];
2590:       mat_i    = c->i[i];
2591:       mat_j    = PetscSafePointerPlusOffset(c->j, mat_i);
2592:       mat_a    = PetscSafePointerPlusOffset(c_a, mat_i);
2593:       mat_ilen = c->ilen + i;
2594:       for (k = kstart; k < kend; k++) {
2595:         if ((tcol = smap[a->j[k]])) {
2596:           *mat_j++ = tcol - 1;
2597:           *mat_a++ = aa[k];
2598:           (*mat_ilen)++;
2599:         }
2600:       }
2601:     }
2602:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2603:     /* Free work space */
2604:     PetscCall(ISRestoreIndices(iscol, &icol));
2605:     PetscCall(PetscFree(smap));
2606:     PetscCall(PetscFree(lens));
2607:     /* sort */
2608:     for (i = 0; i < nrows; i++) {
2609:       PetscInt ilen;

2611:       mat_i = c->i[i];
2612:       mat_j = PetscSafePointerPlusOffset(c->j, mat_i);
2613:       mat_a = PetscSafePointerPlusOffset(c_a, mat_i);
2614:       ilen  = c->ilen[i];
2615:       PetscCall(PetscSortIntWithScalarArray(ilen, mat_j, mat_a));
2616:     }
2617:     PetscCall(MatSeqAIJRestoreArrayWrite(C, &c_a));
2618:   }
2619: #if defined(PETSC_HAVE_DEVICE)
2620:   PetscCall(MatBindToCPU(C, A->boundtocpu));
2621: #endif
2622:   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2623:   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));

2625:   PetscCall(ISRestoreIndices(isrow, &irow));
2626:   *B = C;
2627:   PetscFunctionReturn(PETSC_SUCCESS);
2628: }

2630: static PetscErrorCode MatGetMultiProcBlock_SeqAIJ(Mat mat, MPI_Comm subComm, MatReuse scall, Mat *subMat)
2631: {
2632:   Mat B;

2634:   PetscFunctionBegin;
2635:   if (scall == MAT_INITIAL_MATRIX) {
2636:     PetscCall(MatCreate(subComm, &B));
2637:     PetscCall(MatSetSizes(B, mat->rmap->n, mat->cmap->n, mat->rmap->n, mat->cmap->n));
2638:     PetscCall(MatSetBlockSizesFromMats(B, mat, mat));
2639:     PetscCall(MatSetType(B, MATSEQAIJ));
2640:     PetscCall(MatDuplicateNoCreate_SeqAIJ(B, mat, MAT_COPY_VALUES, PETSC_TRUE));
2641:     *subMat = B;
2642:   } else {
2643:     PetscCall(MatCopy_SeqAIJ(mat, *subMat, SAME_NONZERO_PATTERN));
2644:   }
2645:   PetscFunctionReturn(PETSC_SUCCESS);
2646: }

2648: static PetscErrorCode MatILUFactor_SeqAIJ(Mat inA, IS row, IS col, const MatFactorInfo *info)
2649: {
2650:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)inA->data;
2651:   Mat         outA;
2652:   PetscBool   row_identity, col_identity;

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

2657:   PetscCall(ISIdentity(row, &row_identity));
2658:   PetscCall(ISIdentity(col, &col_identity));

2660:   outA             = inA;
2661:   outA->factortype = MAT_FACTOR_LU;
2662:   PetscCall(PetscFree(inA->solvertype));
2663:   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &inA->solvertype));

2665:   PetscCall(PetscObjectReference((PetscObject)row));
2666:   PetscCall(ISDestroy(&a->row));

2668:   a->row = row;

2670:   PetscCall(PetscObjectReference((PetscObject)col));
2671:   PetscCall(ISDestroy(&a->col));

2673:   a->col = col;

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

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

2683:   PetscCall(MatMarkDiagonal_SeqAIJ(inA));
2684:   if (row_identity && col_identity) {
2685:     PetscCall(MatLUFactorNumeric_SeqAIJ_inplace(outA, inA, info));
2686:   } else {
2687:     PetscCall(MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(outA, inA, info));
2688:   }
2689:   PetscFunctionReturn(PETSC_SUCCESS);
2690: }

2692: PetscErrorCode MatScale_SeqAIJ(Mat inA, PetscScalar alpha)
2693: {
2694:   Mat_SeqAIJ  *a = (Mat_SeqAIJ *)inA->data;
2695:   PetscScalar *v;
2696:   PetscBLASInt one = 1, bnz;

2698:   PetscFunctionBegin;
2699:   PetscCall(MatSeqAIJGetArray(inA, &v));
2700:   PetscCall(PetscBLASIntCast(a->nz, &bnz));
2701:   PetscCallBLAS("BLASscal", BLASscal_(&bnz, &alpha, v, &one));
2702:   PetscCall(PetscLogFlops(a->nz));
2703:   PetscCall(MatSeqAIJRestoreArray(inA, &v));
2704:   PetscCall(MatSeqAIJInvalidateDiagonal(inA));
2705:   PetscFunctionReturn(PETSC_SUCCESS);
2706: }

2708: PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *submatj)
2709: {
2710:   PetscInt i;

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

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

2719:     if (submatj->rbuf1) {
2720:       PetscCall(PetscFree(submatj->rbuf1[0]));
2721:       PetscCall(PetscFree(submatj->rbuf1));
2722:     }

2724:     for (i = 0; i < submatj->nrqs; ++i) PetscCall(PetscFree(submatj->rbuf3[i]));
2725:     PetscCall(PetscFree3(submatj->req_source2, submatj->rbuf2, submatj->rbuf3));
2726:     PetscCall(PetscFree(submatj->pa));
2727:   }

2729: #if defined(PETSC_USE_CTABLE)
2730:   PetscCall(PetscHMapIDestroy(&submatj->rmap));
2731:   if (submatj->cmap_loc) PetscCall(PetscFree(submatj->cmap_loc));
2732:   PetscCall(PetscFree(submatj->rmap_loc));
2733: #else
2734:   PetscCall(PetscFree(submatj->rmap));
2735: #endif

2737:   if (!submatj->allcolumns) {
2738: #if defined(PETSC_USE_CTABLE)
2739:     PetscCall(PetscHMapIDestroy((PetscHMapI *)&submatj->cmap));
2740: #else
2741:     PetscCall(PetscFree(submatj->cmap));
2742: #endif
2743:   }
2744:   PetscCall(PetscFree(submatj->row2proc));

2746:   PetscCall(PetscFree(submatj));
2747:   PetscFunctionReturn(PETSC_SUCCESS);
2748: }

2750: PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat C)
2751: {
2752:   Mat_SeqAIJ  *c       = (Mat_SeqAIJ *)C->data;
2753:   Mat_SubSppt *submatj = c->submatis1;

2755:   PetscFunctionBegin;
2756:   PetscCall((*submatj->destroy)(C));
2757:   PetscCall(MatDestroySubMatrix_Private(submatj));
2758:   PetscFunctionReturn(PETSC_SUCCESS);
2759: }

2761: /* Note this has code duplication with MatDestroySubMatrices_SeqBAIJ() */
2762: static PetscErrorCode MatDestroySubMatrices_SeqAIJ(PetscInt n, Mat *mat[])
2763: {
2764:   PetscInt     i;
2765:   Mat          C;
2766:   Mat_SeqAIJ  *c;
2767:   Mat_SubSppt *submatj;

2769:   PetscFunctionBegin;
2770:   for (i = 0; i < n; i++) {
2771:     C       = (*mat)[i];
2772:     c       = (Mat_SeqAIJ *)C->data;
2773:     submatj = c->submatis1;
2774:     if (submatj) {
2775:       if (--((PetscObject)C)->refct <= 0) {
2776:         PetscCall(PetscFree(C->factorprefix));
2777:         PetscCall((*submatj->destroy)(C));
2778:         PetscCall(MatDestroySubMatrix_Private(submatj));
2779:         PetscCall(PetscFree(C->defaultvectype));
2780:         PetscCall(PetscFree(C->defaultrandtype));
2781:         PetscCall(PetscLayoutDestroy(&C->rmap));
2782:         PetscCall(PetscLayoutDestroy(&C->cmap));
2783:         PetscCall(PetscHeaderDestroy(&C));
2784:       }
2785:     } else {
2786:       PetscCall(MatDestroy(&C));
2787:     }
2788:   }

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

2793:   PetscCall(PetscFree(*mat));
2794:   PetscFunctionReturn(PETSC_SUCCESS);
2795: }

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

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

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

2808: static PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A, PetscInt is_max, IS is[], PetscInt ov)
2809: {
2810:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
2811:   PetscInt        row, i, j, k, l, ll, m, n, *nidx, isz, val;
2812:   const PetscInt *idx;
2813:   PetscInt        start, end, *ai, *aj, bs = (A->rmap->bs > 0 && A->rmap->bs == A->cmap->bs) ? A->rmap->bs : 1;
2814:   PetscBT         table;

2816:   PetscFunctionBegin;
2817:   m  = A->rmap->n / bs;
2818:   ai = a->i;
2819:   aj = a->j;

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

2823:   PetscCall(PetscMalloc1(m + 1, &nidx));
2824:   PetscCall(PetscBTCreate(m, &table));

2826:   for (i = 0; i < is_max; i++) {
2827:     /* Initialize the two local arrays */
2828:     isz = 0;
2829:     PetscCall(PetscBTMemzero(m, table));

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

2835:     if (bs > 1) {
2836:       /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
2837:       for (j = 0; j < n; ++j) {
2838:         if (!PetscBTLookupSet(table, idx[j] / bs)) nidx[isz++] = idx[j] / bs;
2839:       }
2840:       PetscCall(ISRestoreIndices(is[i], &idx));
2841:       PetscCall(ISDestroy(&is[i]));

2843:       k = 0;
2844:       for (j = 0; j < ov; j++) { /* for each overlap */
2845:         n = isz;
2846:         for (; k < n; k++) { /* do only those rows in nidx[k], which are not done yet */
2847:           for (ll = 0; ll < bs; ll++) {
2848:             row   = bs * nidx[k] + ll;
2849:             start = ai[row];
2850:             end   = ai[row + 1];
2851:             for (l = start; l < end; l++) {
2852:               val = aj[l] / bs;
2853:               if (!PetscBTLookupSet(table, val)) nidx[isz++] = val;
2854:             }
2855:           }
2856:         }
2857:       }
2858:       PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, isz, nidx, PETSC_COPY_VALUES, is + i));
2859:     } else {
2860:       /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
2861:       for (j = 0; j < n; ++j) {
2862:         if (!PetscBTLookupSet(table, idx[j])) nidx[isz++] = idx[j];
2863:       }
2864:       PetscCall(ISRestoreIndices(is[i], &idx));
2865:       PetscCall(ISDestroy(&is[i]));

2867:       k = 0;
2868:       for (j = 0; j < ov; j++) { /* for each overlap */
2869:         n = isz;
2870:         for (; k < n; k++) { /* do only those rows in nidx[k], which are not done yet */
2871:           row   = nidx[k];
2872:           start = ai[row];
2873:           end   = ai[row + 1];
2874:           for (l = start; l < end; l++) {
2875:             val = aj[l];
2876:             if (!PetscBTLookupSet(table, val)) nidx[isz++] = val;
2877:           }
2878:         }
2879:       }
2880:       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, PETSC_COPY_VALUES, is + i));
2881:     }
2882:   }
2883:   PetscCall(PetscBTDestroy(&table));
2884:   PetscCall(PetscFree(nidx));
2885:   PetscFunctionReturn(PETSC_SUCCESS);
2886: }

2888: static PetscErrorCode MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
2889: {
2890:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
2891:   PetscInt        i, nz = 0, m = A->rmap->n, n = A->cmap->n;
2892:   const PetscInt *row, *col;
2893:   PetscInt       *cnew, j, *lens;
2894:   IS              icolp, irowp;
2895:   PetscInt       *cwork = NULL;
2896:   PetscScalar    *vwork = NULL;

2898:   PetscFunctionBegin;
2899:   PetscCall(ISInvertPermutation(rowp, PETSC_DECIDE, &irowp));
2900:   PetscCall(ISGetIndices(irowp, &row));
2901:   PetscCall(ISInvertPermutation(colp, PETSC_DECIDE, &icolp));
2902:   PetscCall(ISGetIndices(icolp, &col));

2904:   /* determine lengths of permuted rows */
2905:   PetscCall(PetscMalloc1(m + 1, &lens));
2906:   for (i = 0; i < m; i++) lens[row[i]] = a->i[i + 1] - a->i[i];
2907:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
2908:   PetscCall(MatSetSizes(*B, m, n, m, n));
2909:   PetscCall(MatSetBlockSizesFromMats(*B, A, A));
2910:   PetscCall(MatSetType(*B, ((PetscObject)A)->type_name));
2911:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*B, 0, lens));
2912:   PetscCall(PetscFree(lens));

2914:   PetscCall(PetscMalloc1(n, &cnew));
2915:   for (i = 0; i < m; i++) {
2916:     PetscCall(MatGetRow_SeqAIJ(A, i, &nz, &cwork, &vwork));
2917:     for (j = 0; j < nz; j++) cnew[j] = col[cwork[j]];
2918:     PetscCall(MatSetValues_SeqAIJ(*B, 1, &row[i], nz, cnew, vwork, INSERT_VALUES));
2919:     PetscCall(MatRestoreRow_SeqAIJ(A, i, &nz, &cwork, &vwork));
2920:   }
2921:   PetscCall(PetscFree(cnew));

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

2925: #if defined(PETSC_HAVE_DEVICE)
2926:   PetscCall(MatBindToCPU(*B, A->boundtocpu));
2927: #endif
2928:   PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY));
2929:   PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY));
2930:   PetscCall(ISRestoreIndices(irowp, &row));
2931:   PetscCall(ISRestoreIndices(icolp, &col));
2932:   PetscCall(ISDestroy(&irowp));
2933:   PetscCall(ISDestroy(&icolp));
2934:   if (rowp == colp) PetscCall(MatPropagateSymmetryOptions(A, *B));
2935:   PetscFunctionReturn(PETSC_SUCCESS);
2936: }

2938: PetscErrorCode MatCopy_SeqAIJ(Mat A, Mat B, MatStructure str)
2939: {
2940:   PetscFunctionBegin;
2941:   /* If the two matrices have the same copy implementation, use fast copy. */
2942:   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
2943:     Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
2944:     Mat_SeqAIJ        *b = (Mat_SeqAIJ *)B->data;
2945:     const PetscScalar *aa;

2947:     PetscCall(MatSeqAIJGetArrayRead(A, &aa));
2948:     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]);
2949:     PetscCall(PetscArraycpy(b->a, aa, a->i[A->rmap->n]));
2950:     PetscCall(PetscObjectStateIncrease((PetscObject)B));
2951:     PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
2952:   } else {
2953:     PetscCall(MatCopy_Basic(A, B, str));
2954:   }
2955:   PetscFunctionReturn(PETSC_SUCCESS);
2956: }

2958: PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat A, PetscScalar *array[])
2959: {
2960:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

2962:   PetscFunctionBegin;
2963:   *array = a->a;
2964:   PetscFunctionReturn(PETSC_SUCCESS);
2965: }

2967: PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat A, PetscScalar *array[])
2968: {
2969:   PetscFunctionBegin;
2970:   *array = NULL;
2971:   PetscFunctionReturn(PETSC_SUCCESS);
2972: }

2974: /*
2975:    Computes the number of nonzeros per row needed for preallocation when X and Y
2976:    have different nonzero structure.
2977: */
2978: PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt m, const PetscInt *xi, const PetscInt *xj, const PetscInt *yi, const PetscInt *yj, PetscInt *nnz)
2979: {
2980:   PetscInt i, j, k, nzx, nzy;

2982:   PetscFunctionBegin;
2983:   /* Set the number of nonzeros in the new matrix */
2984:   for (i = 0; i < m; i++) {
2985:     const PetscInt *xjj = PetscSafePointerPlusOffset(xj, xi[i]), *yjj = PetscSafePointerPlusOffset(yj, yi[i]);
2986:     nzx    = xi[i + 1] - xi[i];
2987:     nzy    = yi[i + 1] - yi[i];
2988:     nnz[i] = 0;
2989:     for (j = 0, k = 0; j < nzx; j++) {                  /* Point in X */
2990:       for (; k < nzy && yjj[k] < xjj[j]; k++) nnz[i]++; /* Catch up to X */
2991:       if (k < nzy && yjj[k] == xjj[j]) k++;             /* Skip duplicate */
2992:       nnz[i]++;
2993:     }
2994:     for (; k < nzy; k++) nnz[i]++;
2995:   }
2996:   PetscFunctionReturn(PETSC_SUCCESS);
2997: }

2999: PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat Y, Mat X, PetscInt *nnz)
3000: {
3001:   PetscInt    m = Y->rmap->N;
3002:   Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data;
3003:   Mat_SeqAIJ *y = (Mat_SeqAIJ *)Y->data;

3005:   PetscFunctionBegin;
3006:   /* Set the number of nonzeros in the new matrix */
3007:   PetscCall(MatAXPYGetPreallocation_SeqX_private(m, x->i, x->j, y->i, y->j, nnz));
3008:   PetscFunctionReturn(PETSC_SUCCESS);
3009: }

3011: PetscErrorCode MatAXPY_SeqAIJ(Mat Y, PetscScalar a, Mat X, MatStructure str)
3012: {
3013:   Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data, *y = (Mat_SeqAIJ *)Y->data;

3015:   PetscFunctionBegin;
3016:   if (str == UNKNOWN_NONZERO_PATTERN || (PetscDefined(USE_DEBUG) && str == SAME_NONZERO_PATTERN)) {
3017:     PetscBool e = x->nz == y->nz ? PETSC_TRUE : PETSC_FALSE;
3018:     if (e) {
3019:       PetscCall(PetscArraycmp(x->i, y->i, Y->rmap->n + 1, &e));
3020:       if (e) {
3021:         PetscCall(PetscArraycmp(x->j, y->j, y->nz, &e));
3022:         if (e) str = SAME_NONZERO_PATTERN;
3023:       }
3024:     }
3025:     if (!e) PetscCheck(str != SAME_NONZERO_PATTERN, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatStructure is not SAME_NONZERO_PATTERN");
3026:   }
3027:   if (str == SAME_NONZERO_PATTERN) {
3028:     const PetscScalar *xa;
3029:     PetscScalar       *ya, alpha = a;
3030:     PetscBLASInt       one = 1, bnz;

3032:     PetscCall(PetscBLASIntCast(x->nz, &bnz));
3033:     PetscCall(MatSeqAIJGetArray(Y, &ya));
3034:     PetscCall(MatSeqAIJGetArrayRead(X, &xa));
3035:     PetscCallBLAS("BLASaxpy", BLASaxpy_(&bnz, &alpha, xa, &one, ya, &one));
3036:     PetscCall(MatSeqAIJRestoreArrayRead(X, &xa));
3037:     PetscCall(MatSeqAIJRestoreArray(Y, &ya));
3038:     PetscCall(PetscLogFlops(2.0 * bnz));
3039:     PetscCall(MatSeqAIJInvalidateDiagonal(Y));
3040:     PetscCall(PetscObjectStateIncrease((PetscObject)Y));
3041:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
3042:     PetscCall(MatAXPY_Basic(Y, a, X, str));
3043:   } else {
3044:     Mat       B;
3045:     PetscInt *nnz;
3046:     PetscCall(PetscMalloc1(Y->rmap->N, &nnz));
3047:     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &B));
3048:     PetscCall(PetscObjectSetName((PetscObject)B, ((PetscObject)Y)->name));
3049:     PetscCall(MatSetLayouts(B, Y->rmap, Y->cmap));
3050:     PetscCall(MatSetType(B, ((PetscObject)Y)->type_name));
3051:     PetscCall(MatAXPYGetPreallocation_SeqAIJ(Y, X, nnz));
3052:     PetscCall(MatSeqAIJSetPreallocation(B, 0, nnz));
3053:     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
3054:     PetscCall(MatHeaderMerge(Y, &B));
3055:     PetscCall(MatSeqAIJCheckInode(Y));
3056:     PetscCall(PetscFree(nnz));
3057:   }
3058:   PetscFunctionReturn(PETSC_SUCCESS);
3059: }

3061: PETSC_INTERN PetscErrorCode MatConjugate_SeqAIJ(Mat mat)
3062: {
3063: #if defined(PETSC_USE_COMPLEX)
3064:   Mat_SeqAIJ  *aij = (Mat_SeqAIJ *)mat->data;
3065:   PetscInt     i, nz;
3066:   PetscScalar *a;

3068:   PetscFunctionBegin;
3069:   nz = aij->nz;
3070:   PetscCall(MatSeqAIJGetArray(mat, &a));
3071:   for (i = 0; i < nz; i++) a[i] = PetscConj(a[i]);
3072:   PetscCall(MatSeqAIJRestoreArray(mat, &a));
3073: #else
3074:   PetscFunctionBegin;
3075: #endif
3076:   PetscFunctionReturn(PETSC_SUCCESS);
3077: }

3079: static PetscErrorCode MatGetRowMaxAbs_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3080: {
3081:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3082:   PetscInt         i, j, m = A->rmap->n, *ai, *aj, ncols, n;
3083:   PetscReal        atmp;
3084:   PetscScalar     *x;
3085:   const MatScalar *aa, *av;

3087:   PetscFunctionBegin;
3088:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3089:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3090:   aa = av;
3091:   ai = a->i;
3092:   aj = a->j;

3094:   PetscCall(VecGetArrayWrite(v, &x));
3095:   PetscCall(VecGetLocalSize(v, &n));
3096:   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3097:   for (i = 0; i < m; i++) {
3098:     ncols = ai[1] - ai[0];
3099:     ai++;
3100:     x[i] = 0;
3101:     for (j = 0; j < ncols; j++) {
3102:       atmp = PetscAbsScalar(*aa);
3103:       if (PetscAbsScalar(x[i]) < atmp) {
3104:         x[i] = atmp;
3105:         if (idx) idx[i] = *aj;
3106:       }
3107:       aa++;
3108:       aj++;
3109:     }
3110:   }
3111:   PetscCall(VecRestoreArrayWrite(v, &x));
3112:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3113:   PetscFunctionReturn(PETSC_SUCCESS);
3114: }

3116: static PetscErrorCode MatGetRowSumAbs_SeqAIJ(Mat A, Vec v)
3117: {
3118:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3119:   PetscInt         i, j, m = A->rmap->n, *ai, ncols, n;
3120:   PetscScalar     *x;
3121:   const MatScalar *aa, *av;

3123:   PetscFunctionBegin;
3124:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3125:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3126:   aa = av;
3127:   ai = a->i;

3129:   PetscCall(VecGetArrayWrite(v, &x));
3130:   PetscCall(VecGetLocalSize(v, &n));
3131:   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3132:   for (i = 0; i < m; i++) {
3133:     ncols = ai[1] - ai[0];
3134:     ai++;
3135:     x[i] = 0;
3136:     for (j = 0; j < ncols; j++) {
3137:       x[i] += PetscAbsScalar(*aa);
3138:       aa++;
3139:     }
3140:   }
3141:   PetscCall(VecRestoreArrayWrite(v, &x));
3142:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3143:   PetscFunctionReturn(PETSC_SUCCESS);
3144: }

3146: static PetscErrorCode MatGetRowMax_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3147: {
3148:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3149:   PetscInt         i, j, m = A->rmap->n, *ai, *aj, ncols, n;
3150:   PetscScalar     *x;
3151:   const MatScalar *aa, *av;

3153:   PetscFunctionBegin;
3154:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3155:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3156:   aa = av;
3157:   ai = a->i;
3158:   aj = a->j;

3160:   PetscCall(VecGetArrayWrite(v, &x));
3161:   PetscCall(VecGetLocalSize(v, &n));
3162:   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3163:   for (i = 0; i < m; i++) {
3164:     ncols = ai[1] - ai[0];
3165:     ai++;
3166:     if (ncols == A->cmap->n) { /* row is dense */
3167:       x[i] = *aa;
3168:       if (idx) idx[i] = 0;
3169:     } else { /* row is sparse so already KNOW maximum is 0.0 or higher */
3170:       x[i] = 0.0;
3171:       if (idx) {
3172:         for (j = 0; j < ncols; j++) { /* find first implicit 0.0 in the row */
3173:           if (aj[j] > j) {
3174:             idx[i] = j;
3175:             break;
3176:           }
3177:         }
3178:         /* in case first implicit 0.0 in the row occurs at ncols-th column */
3179:         if (j == ncols && j < A->cmap->n) idx[i] = j;
3180:       }
3181:     }
3182:     for (j = 0; j < ncols; j++) {
3183:       if (PetscRealPart(x[i]) < PetscRealPart(*aa)) {
3184:         x[i] = *aa;
3185:         if (idx) idx[i] = *aj;
3186:       }
3187:       aa++;
3188:       aj++;
3189:     }
3190:   }
3191:   PetscCall(VecRestoreArrayWrite(v, &x));
3192:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3193:   PetscFunctionReturn(PETSC_SUCCESS);
3194: }

3196: static PetscErrorCode MatGetRowMinAbs_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3197: {
3198:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3199:   PetscInt         i, j, m = A->rmap->n, *ai, *aj, ncols, n;
3200:   PetscScalar     *x;
3201:   const MatScalar *aa, *av;

3203:   PetscFunctionBegin;
3204:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3205:   aa = av;
3206:   ai = a->i;
3207:   aj = a->j;

3209:   PetscCall(VecGetArrayWrite(v, &x));
3210:   PetscCall(VecGetLocalSize(v, &n));
3211:   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector, %" PetscInt_FMT " vs. %" PetscInt_FMT " rows", m, n);
3212:   for (i = 0; i < m; i++) {
3213:     ncols = ai[1] - ai[0];
3214:     ai++;
3215:     if (ncols == A->cmap->n) { /* row is dense */
3216:       x[i] = *aa;
3217:       if (idx) idx[i] = 0;
3218:     } else { /* row is sparse so already KNOW minimum is 0.0 or higher */
3219:       x[i] = 0.0;
3220:       if (idx) { /* find first implicit 0.0 in the row */
3221:         for (j = 0; j < ncols; j++) {
3222:           if (aj[j] > j) {
3223:             idx[i] = j;
3224:             break;
3225:           }
3226:         }
3227:         /* in case first implicit 0.0 in the row occurs at ncols-th column */
3228:         if (j == ncols && j < A->cmap->n) idx[i] = j;
3229:       }
3230:     }
3231:     for (j = 0; j < ncols; j++) {
3232:       if (PetscAbsScalar(x[i]) > PetscAbsScalar(*aa)) {
3233:         x[i] = *aa;
3234:         if (idx) idx[i] = *aj;
3235:       }
3236:       aa++;
3237:       aj++;
3238:     }
3239:   }
3240:   PetscCall(VecRestoreArrayWrite(v, &x));
3241:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3242:   PetscFunctionReturn(PETSC_SUCCESS);
3243: }

3245: static PetscErrorCode MatGetRowMin_SeqAIJ(Mat A, Vec v, PetscInt idx[])
3246: {
3247:   Mat_SeqAIJ      *a = (Mat_SeqAIJ *)A->data;
3248:   PetscInt         i, j, m = A->rmap->n, ncols, n;
3249:   const PetscInt  *ai, *aj;
3250:   PetscScalar     *x;
3251:   const MatScalar *aa, *av;

3253:   PetscFunctionBegin;
3254:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3255:   PetscCall(MatSeqAIJGetArrayRead(A, &av));
3256:   aa = av;
3257:   ai = a->i;
3258:   aj = a->j;

3260:   PetscCall(VecGetArrayWrite(v, &x));
3261:   PetscCall(VecGetLocalSize(v, &n));
3262:   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
3263:   for (i = 0; i < m; i++) {
3264:     ncols = ai[1] - ai[0];
3265:     ai++;
3266:     if (ncols == A->cmap->n) { /* row is dense */
3267:       x[i] = *aa;
3268:       if (idx) idx[i] = 0;
3269:     } else { /* row is sparse so already KNOW minimum is 0.0 or lower */
3270:       x[i] = 0.0;
3271:       if (idx) { /* find first implicit 0.0 in the row */
3272:         for (j = 0; j < ncols; j++) {
3273:           if (aj[j] > j) {
3274:             idx[i] = j;
3275:             break;
3276:           }
3277:         }
3278:         /* in case first implicit 0.0 in the row occurs at ncols-th column */
3279:         if (j == ncols && j < A->cmap->n) idx[i] = j;
3280:       }
3281:     }
3282:     for (j = 0; j < ncols; j++) {
3283:       if (PetscRealPart(x[i]) > PetscRealPart(*aa)) {
3284:         x[i] = *aa;
3285:         if (idx) idx[i] = *aj;
3286:       }
3287:       aa++;
3288:       aj++;
3289:     }
3290:   }
3291:   PetscCall(VecRestoreArrayWrite(v, &x));
3292:   PetscCall(MatSeqAIJRestoreArrayRead(A, &av));
3293:   PetscFunctionReturn(PETSC_SUCCESS);
3294: }

3296: static PetscErrorCode MatInvertBlockDiagonal_SeqAIJ(Mat A, const PetscScalar **values)
3297: {
3298:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
3299:   PetscInt        i, bs = PetscAbs(A->rmap->bs), mbs = A->rmap->n / bs, ipvt[5], bs2 = bs * bs, *v_pivots, ij[7], *IJ, j;
3300:   MatScalar      *diag, work[25], *v_work;
3301:   const PetscReal shift = 0.0;
3302:   PetscBool       allowzeropivot, zeropivotdetected = PETSC_FALSE;

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

3427: static PetscErrorCode MatSetRandom_SeqAIJ(Mat x, PetscRandom rctx)
3428: {
3429:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)x->data;
3430:   PetscScalar a, *aa;
3431:   PetscInt    m, n, i, j, col;

3433:   PetscFunctionBegin;
3434:   if (!x->assembled) {
3435:     PetscCall(MatGetSize(x, &m, &n));
3436:     for (i = 0; i < m; i++) {
3437:       for (j = 0; j < aij->imax[i]; j++) {
3438:         PetscCall(PetscRandomGetValue(rctx, &a));
3439:         col = (PetscInt)(n * PetscRealPart(a));
3440:         PetscCall(MatSetValues(x, 1, &i, 1, &col, &a, ADD_VALUES));
3441:       }
3442:     }
3443:   } else {
3444:     PetscCall(MatSeqAIJGetArrayWrite(x, &aa));
3445:     for (i = 0; i < aij->nz; i++) PetscCall(PetscRandomGetValue(rctx, aa + i));
3446:     PetscCall(MatSeqAIJRestoreArrayWrite(x, &aa));
3447:   }
3448:   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
3449:   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
3450:   PetscFunctionReturn(PETSC_SUCCESS);
3451: }

3453: /* Like MatSetRandom_SeqAIJ, but do not set values on columns in range of [low, high) */
3454: PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat x, PetscInt low, PetscInt high, PetscRandom rctx)
3455: {
3456:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)x->data;
3457:   PetscScalar a;
3458:   PetscInt    m, n, i, j, col, nskip;

3460:   PetscFunctionBegin;
3461:   nskip = high - low;
3462:   PetscCall(MatGetSize(x, &m, &n));
3463:   n -= nskip; /* shrink number of columns where nonzeros can be set */
3464:   for (i = 0; i < m; i++) {
3465:     for (j = 0; j < aij->imax[i]; j++) {
3466:       PetscCall(PetscRandomGetValue(rctx, &a));
3467:       col = (PetscInt)(n * PetscRealPart(a));
3468:       if (col >= low) col += nskip; /* shift col rightward to skip the hole */
3469:       PetscCall(MatSetValues(x, 1, &i, 1, &col, &a, ADD_VALUES));
3470:     }
3471:   }
3472:   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
3473:   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
3474:   PetscFunctionReturn(PETSC_SUCCESS);
3475: }

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

3634: static PetscErrorCode MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat, PetscInt *indices)
3635: {
3636:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
3637:   PetscInt    i, nz, n;

3639:   PetscFunctionBegin;
3640:   nz = aij->maxnz;
3641:   n  = mat->rmap->n;
3642:   for (i = 0; i < nz; i++) aij->j[i] = indices[i];
3643:   aij->nz = nz;
3644:   for (i = 0; i < n; i++) aij->ilen[i] = aij->imax[i];
3645:   PetscFunctionReturn(PETSC_SUCCESS);
3646: }

3648: /*
3649:  * Given a sparse matrix with global column indices, compact it by using a local column space.
3650:  * The result matrix helps saving memory in other algorithms, such as MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable()
3651:  */
3652: PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat mat, ISLocalToGlobalMapping *mapping)
3653: {
3654:   Mat_SeqAIJ   *aij = (Mat_SeqAIJ *)mat->data;
3655:   PetscHMapI    gid1_lid1;
3656:   PetscHashIter tpos;
3657:   PetscInt      gid, lid, i, ec, nz = aij->nz;
3658:   PetscInt     *garray, *jj = aij->j;

3660:   PetscFunctionBegin;
3662:   PetscAssertPointer(mapping, 2);
3663:   /* use a table */
3664:   PetscCall(PetscHMapICreateWithSize(mat->rmap->n, &gid1_lid1));
3665:   ec = 0;
3666:   for (i = 0; i < nz; i++) {
3667:     PetscInt data, gid1 = jj[i] + 1;
3668:     PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &data));
3669:     if (!data) {
3670:       /* one based table */
3671:       PetscCall(PetscHMapISet(gid1_lid1, gid1, ++ec));
3672:     }
3673:   }
3674:   /* form array of columns we need */
3675:   PetscCall(PetscMalloc1(ec, &garray));
3676:   PetscHashIterBegin(gid1_lid1, tpos);
3677:   while (!PetscHashIterAtEnd(gid1_lid1, tpos)) {
3678:     PetscHashIterGetKey(gid1_lid1, tpos, gid);
3679:     PetscHashIterGetVal(gid1_lid1, tpos, lid);
3680:     PetscHashIterNext(gid1_lid1, tpos);
3681:     gid--;
3682:     lid--;
3683:     garray[lid] = gid;
3684:   }
3685:   PetscCall(PetscSortInt(ec, garray)); /* sort, and rebuild */
3686:   PetscCall(PetscHMapIClear(gid1_lid1));
3687:   for (i = 0; i < ec; i++) PetscCall(PetscHMapISet(gid1_lid1, garray[i] + 1, i + 1));
3688:   /* compact out the extra columns in B */
3689:   for (i = 0; i < nz; i++) {
3690:     PetscInt gid1 = jj[i] + 1;
3691:     PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &lid));
3692:     lid--;
3693:     jj[i] = lid;
3694:   }
3695:   PetscCall(PetscLayoutDestroy(&mat->cmap));
3696:   PetscCall(PetscHMapIDestroy(&gid1_lid1));
3697:   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)mat), ec, ec, 1, &mat->cmap));
3698:   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, mat->cmap->bs, mat->cmap->n, garray, PETSC_OWN_POINTER, mapping));
3699:   PetscCall(ISLocalToGlobalMappingSetType(*mapping, ISLOCALTOGLOBALMAPPINGHASH));
3700:   PetscFunctionReturn(PETSC_SUCCESS);
3701: }

3703: /*@
3704:   MatSeqAIJSetColumnIndices - Set the column indices for all the rows
3705:   in the matrix.

3707:   Input Parameters:
3708: + mat     - the `MATSEQAIJ` matrix
3709: - indices - the column indices

3711:   Level: advanced

3713:   Notes:
3714:   This can be called if you have precomputed the nonzero structure of the
3715:   matrix and want to provide it to the matrix object to improve the performance
3716:   of the `MatSetValues()` operation.

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

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

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

3725: .seealso: [](ch_matrices), `Mat`, `MATSEQAIJ`
3726: @*/
3727: PetscErrorCode MatSeqAIJSetColumnIndices(Mat mat, PetscInt *indices)
3728: {
3729:   PetscFunctionBegin;
3731:   PetscAssertPointer(indices, 2);
3732:   PetscUseMethod(mat, "MatSeqAIJSetColumnIndices_C", (Mat, PetscInt *), (mat, indices));
3733:   PetscFunctionReturn(PETSC_SUCCESS);
3734: }

3736: static PetscErrorCode MatStoreValues_SeqAIJ(Mat mat)
3737: {
3738:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
3739:   size_t      nz  = aij->i[mat->rmap->n];

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

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

3747:   /* copy values over */
3748:   PetscCall(PetscArraycpy(aij->saved_values, aij->a, nz));
3749:   PetscFunctionReturn(PETSC_SUCCESS);
3750: }

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

3756:   Logically Collect

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

3761:   Level: advanced

3763:   Example Usage:
3764: .vb
3765:     Using SNES
3766:     Create Jacobian matrix
3767:     Set linear terms into matrix
3768:     Apply boundary conditions to matrix, at this time matrix must have
3769:       final nonzero structure (i.e. setting the nonlinear terms and applying
3770:       boundary conditions again will not change the nonzero structure
3771:     MatSetOption(mat, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE);
3772:     MatStoreValues(mat);
3773:     Call SNESSetJacobian() with matrix
3774:     In your Jacobian routine
3775:       MatRetrieveValues(mat);
3776:       Set nonlinear terms in matrix

3778:     Without `SNESSolve()`, i.e. when you handle nonlinear solve yourself:
3779:     // build linear portion of Jacobian
3780:     MatSetOption(mat, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE);
3781:     MatStoreValues(mat);
3782:     loop over nonlinear iterations
3783:        MatRetrieveValues(mat);
3784:        // call MatSetValues(mat,...) to set nonliner portion of Jacobian
3785:        // call MatAssemblyBegin/End() on matrix
3786:        Solve linear system with Jacobian
3787:     endloop
3788: .ve

3790:   Notes:
3791:   Matrix must already be assembled before calling this routine
3792:   Must set the matrix option `MatSetOption`(mat,`MAT_NEW_NONZERO_LOCATIONS`,`PETSC_FALSE`); before
3793:   calling this routine.

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

3798: .seealso: [](ch_matrices), `Mat`, `MatRetrieveValues()`
3799: @*/
3800: PetscErrorCode MatStoreValues(Mat mat)
3801: {
3802:   PetscFunctionBegin;
3804:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3805:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3806:   PetscUseMethod(mat, "MatStoreValues_C", (Mat), (mat));
3807:   PetscFunctionReturn(PETSC_SUCCESS);
3808: }

3810: static PetscErrorCode MatRetrieveValues_SeqAIJ(Mat mat)
3811: {
3812:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
3813:   PetscInt    nz  = aij->i[mat->rmap->n];

3815:   PetscFunctionBegin;
3816:   PetscCheck(aij->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
3817:   PetscCheck(aij->saved_values, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatStoreValues(A);first");
3818:   /* copy values over */
3819:   PetscCall(PetscArraycpy(aij->a, aij->saved_values, nz));
3820:   PetscFunctionReturn(PETSC_SUCCESS);
3821: }

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

3826:   Logically Collect

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

3831:   Level: advanced

3833: .seealso: [](ch_matrices), `Mat`, `MatStoreValues()`
3834: @*/
3835: PetscErrorCode MatRetrieveValues(Mat mat)
3836: {
3837:   PetscFunctionBegin;
3839:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3840:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3841:   PetscUseMethod(mat, "MatRetrieveValues_C", (Mat), (mat));
3842:   PetscFunctionReturn(PETSC_SUCCESS);
3843: }

3845: /*@
3846:   MatCreateSeqAIJ - Creates a sparse matrix in `MATSEQAIJ` (compressed row) format
3847:   (the default parallel PETSc format).  For good matrix assembly performance
3848:   the user should preallocate the matrix storage by setting the parameter `nz`
3849:   (or the array `nnz`).

3851:   Collective

3853:   Input Parameters:
3854: + comm - MPI communicator, set to `PETSC_COMM_SELF`
3855: . m    - number of rows
3856: . n    - number of columns
3857: . nz   - number of nonzeros per row (same for all rows)
3858: - nnz  - array containing the number of nonzeros in the various rows
3859:          (possibly different for each row) or NULL

3861:   Output Parameter:
3862: . A - the matrix

3864:   Options Database Keys:
3865: + -mat_no_inode            - Do not use inodes
3866: - -mat_inode_limit <limit> - Sets inode limit (max limit=5)

3868:   Level: intermediate

3870:   Notes:
3871:   It is recommend to use `MatCreateFromOptions()` instead of this routine

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

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

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

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

3889: .seealso: [](ch_matrices), `Mat`, [Sparse Matrix Creation](sec_matsparse), `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()`
3890: @*/
3891: PetscErrorCode MatCreateSeqAIJ(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A)
3892: {
3893:   PetscFunctionBegin;
3894:   PetscCall(MatCreate(comm, A));
3895:   PetscCall(MatSetSizes(*A, m, n, m, n));
3896:   PetscCall(MatSetType(*A, MATSEQAIJ));
3897:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*A, nz, nnz));
3898:   PetscFunctionReturn(PETSC_SUCCESS);
3899: }

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

3907:   Collective

3909:   Input Parameters:
3910: + B   - The matrix
3911: . nz  - number of nonzeros per row (same for all rows)
3912: - nnz - array containing the number of nonzeros in the various rows
3913:          (possibly different for each row) or NULL

3915:   Options Database Keys:
3916: + -mat_no_inode            - Do not use inodes
3917: - -mat_inode_limit <limit> - Sets inode limit (max limit=5)

3919:   Level: intermediate

3921:   Notes:
3922:   If `nnz` is given then `nz` is ignored

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

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

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

3938:   Developer Notes:
3939:   Use nz of `MAT_SKIP_ALLOCATION` to not allocate any space for the matrix
3940:   entries or columns indices

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

3947: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()`, `MatGetInfo()`,
3948:           `MatSeqAIJSetTotalPreallocation()`
3949: @*/
3950: PetscErrorCode MatSeqAIJSetPreallocation(Mat B, PetscInt nz, const PetscInt nnz[])
3951: {
3952:   PetscFunctionBegin;
3955:   PetscTryMethod(B, "MatSeqAIJSetPreallocation_C", (Mat, PetscInt, const PetscInt[]), (B, nz, nnz));
3956:   PetscFunctionReturn(PETSC_SUCCESS);
3957: }

3959: PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat B, PetscInt nz, const PetscInt *nnz)
3960: {
3961:   Mat_SeqAIJ *b              = (Mat_SeqAIJ *)B->data;
3962:   PetscBool   skipallocation = PETSC_FALSE, realalloc = PETSC_FALSE;
3963:   PetscInt    i;

3965:   PetscFunctionBegin;
3966:   if (B->hash_active) {
3967:     B->ops[0] = b->cops;
3968:     PetscCall(PetscHMapIJVDestroy(&b->ht));
3969:     PetscCall(PetscFree(b->dnz));
3970:     B->hash_active = PETSC_FALSE;
3971:   }
3972:   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
3973:   if (nz == MAT_SKIP_ALLOCATION) {
3974:     skipallocation = PETSC_TRUE;
3975:     nz             = 0;
3976:   }
3977:   PetscCall(PetscLayoutSetUp(B->rmap));
3978:   PetscCall(PetscLayoutSetUp(B->cmap));

3980:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
3981:   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nz cannot be less than 0: value %" PetscInt_FMT, nz);
3982:   if (nnz) {
3983:     for (i = 0; i < B->rmap->n; i++) {
3984:       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]);
3985:       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);
3986:     }
3987:   }

3989:   B->preallocated = PETSC_TRUE;
3990:   if (!skipallocation) {
3991:     if (!b->imax) { PetscCall(PetscMalloc1(B->rmap->n, &b->imax)); }
3992:     if (!b->ilen) {
3993:       /* b->ilen will count nonzeros in each row so far. */
3994:       PetscCall(PetscCalloc1(B->rmap->n, &b->ilen));
3995:     } else {
3996:       PetscCall(PetscMemzero(b->ilen, B->rmap->n * sizeof(PetscInt)));
3997:     }
3998:     if (!b->ipre) PetscCall(PetscMalloc1(B->rmap->n, &b->ipre));
3999:     if (!nnz) {
4000:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
4001:       else if (nz < 0) nz = 1;
4002:       nz = PetscMin(nz, B->cmap->n);
4003:       for (i = 0; i < B->rmap->n; i++) b->imax[i] = nz;
4004:       PetscCall(PetscIntMultError(nz, B->rmap->n, &nz));
4005:     } else {
4006:       PetscInt64 nz64 = 0;
4007:       for (i = 0; i < B->rmap->n; i++) {
4008:         b->imax[i] = nnz[i];
4009:         nz64 += nnz[i];
4010:       }
4011:       PetscCall(PetscIntCast(nz64, &nz));
4012:     }

4014:     /* allocate the matrix space */
4015:     PetscCall(MatSeqXAIJFreeAIJ(B, &b->a, &b->j, &b->i));
4016:     PetscCall(PetscShmgetAllocateArray(nz, sizeof(PetscInt), (void **)&b->j));
4017:     PetscCall(PetscShmgetAllocateArray(B->rmap->n + 1, sizeof(PetscInt), (void **)&b->i));
4018:     b->free_ij = PETSC_TRUE;
4019:     if (B->structure_only) {
4020:       b->free_a = PETSC_FALSE;
4021:     } else {
4022:       PetscCall(PetscShmgetAllocateArray(nz, sizeof(PetscScalar), (void **)&b->a));
4023:       b->free_a = PETSC_TRUE;
4024:     }
4025:     b->i[0] = 0;
4026:     for (i = 1; i < B->rmap->n + 1; i++) b->i[i] = b->i[i - 1] + b->imax[i - 1];
4027:   } else {
4028:     b->free_a  = PETSC_FALSE;
4029:     b->free_ij = PETSC_FALSE;
4030:   }

4032:   if (b->ipre && nnz != b->ipre && b->imax) {
4033:     /* reserve user-requested sparsity */
4034:     PetscCall(PetscArraycpy(b->ipre, b->imax, B->rmap->n));
4035:   }

4037:   b->nz               = 0;
4038:   b->maxnz            = nz;
4039:   B->info.nz_unneeded = (double)b->maxnz;
4040:   if (realalloc) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
4041:   B->was_assembled = PETSC_FALSE;
4042:   B->assembled     = PETSC_FALSE;
4043:   /* We simply deem preallocation has changed nonzero state. Updating the state
4044:      will give clients (like AIJKokkos) a chance to know something has happened.
4045:   */
4046:   B->nonzerostate++;
4047:   PetscFunctionReturn(PETSC_SUCCESS);
4048: }

4050: static PetscErrorCode MatResetPreallocation_SeqAIJ(Mat A)
4051: {
4052:   Mat_SeqAIJ *a;
4053:   PetscInt    i;
4054:   PetscBool   skipreset;

4056:   PetscFunctionBegin;

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

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

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

4068:   PetscCall(PetscArraycmp(a->ipre, a->ilen, A->rmap->n, &skipreset));
4069:   if (!skipreset) {
4070:     PetscCall(PetscArraycpy(a->imax, a->ipre, A->rmap->n));
4071:     PetscCall(PetscArrayzero(a->ilen, A->rmap->n));
4072:     a->i[0] = 0;
4073:     for (i = 1; i < A->rmap->n + 1; i++) a->i[i] = a->i[i - 1] + a->imax[i - 1];
4074:     A->preallocated     = PETSC_TRUE;
4075:     a->nz               = 0;
4076:     a->maxnz            = a->i[A->rmap->n];
4077:     A->info.nz_unneeded = (double)a->maxnz;
4078:     A->was_assembled    = PETSC_FALSE;
4079:     A->assembled        = PETSC_FALSE;
4080:     A->nonzerostate++;
4081:     /* Log that the state of this object has changed; this will help guarantee that preconditioners get re-setup */
4082:     PetscCall(PetscObjectStateIncrease((PetscObject)A));
4083:   }
4084:   PetscFunctionReturn(PETSC_SUCCESS);
4085: }

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

4090:   Input Parameters:
4091: + B - the matrix
4092: . i - the indices into `j` for the start of each row (indices start with zero)
4093: . j - the column indices for each row (indices start with zero) these must be sorted for each row
4094: - v - optional values in the matrix, use `NULL` if not provided

4096:   Level: developer

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

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

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

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

4110: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatSeqAIJSetPreallocation()`, `MATSEQAIJ`, `MatResetPreallocation()`
4111: @*/
4112: PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B, const PetscInt i[], const PetscInt j[], const PetscScalar v[])
4113: {
4114:   PetscFunctionBegin;
4117:   PetscTryMethod(B, "MatSeqAIJSetPreallocationCSR_C", (Mat, const PetscInt[], const PetscInt[], const PetscScalar[]), (B, i, j, v));
4118:   PetscFunctionReturn(PETSC_SUCCESS);
4119: }

4121: static PetscErrorCode MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B, const PetscInt Ii[], const PetscInt J[], const PetscScalar v[])
4122: {
4123:   PetscInt  i;
4124:   PetscInt  m, n;
4125:   PetscInt  nz;
4126:   PetscInt *nnz;

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

4131:   PetscCall(PetscLayoutSetUp(B->rmap));
4132:   PetscCall(PetscLayoutSetUp(B->cmap));

4134:   PetscCall(MatGetSize(B, &m, &n));
4135:   PetscCall(PetscMalloc1(m + 1, &nnz));
4136:   for (i = 0; i < m; i++) {
4137:     nz = Ii[i + 1] - Ii[i];
4138:     PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local row %" PetscInt_FMT " has a negative number of columns %" PetscInt_FMT, i, nz);
4139:     nnz[i] = nz;
4140:   }
4141:   PetscCall(MatSeqAIJSetPreallocation(B, 0, nnz));
4142:   PetscCall(PetscFree(nnz));

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

4146:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
4147:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));

4149:   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
4150:   PetscFunctionReturn(PETSC_SUCCESS);
4151: }

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

4156:   Input Parameters:
4157: + A     - left-hand side matrix
4158: . B     - right-hand side matrix
4159: - reuse - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`

4161:   Output Parameter:
4162: . C - Kronecker product of `A` and `B`

4164:   Level: intermediate

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

4169: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MATSEQAIJ`, `MATKAIJ`, `MatReuse`
4170: @*/
4171: PetscErrorCode MatSeqAIJKron(Mat A, Mat B, MatReuse reuse, Mat *C)
4172: {
4173:   PetscFunctionBegin;
4178:   PetscAssertPointer(C, 4);
4179:   if (reuse == MAT_REUSE_MATRIX) {
4182:   }
4183:   PetscTryMethod(A, "MatSeqAIJKron_C", (Mat, Mat, MatReuse, Mat *), (A, B, reuse, C));
4184:   PetscFunctionReturn(PETSC_SUCCESS);
4185: }

4187: static PetscErrorCode MatSeqAIJKron_SeqAIJ(Mat A, Mat B, MatReuse reuse, Mat *C)
4188: {
4189:   Mat                newmat;
4190:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data;
4191:   Mat_SeqAIJ        *b = (Mat_SeqAIJ *)B->data;
4192:   PetscScalar       *v;
4193:   const PetscScalar *aa, *ba;
4194:   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;
4195:   PetscBool          flg;

4197:   PetscFunctionBegin;
4198:   PetscCheck(!A->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4199:   PetscCheck(A->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
4200:   PetscCheck(!B->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4201:   PetscCheck(B->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
4202:   PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &flg));
4203:   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatType %s", ((PetscObject)B)->type_name);
4204:   PetscCheck(reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatReuse %d", (int)reuse);
4205:   if (reuse == MAT_INITIAL_MATRIX) {
4206:     PetscCall(PetscMalloc2(am * bm + 1, &i, a->i[am] * b->i[bm], &j));
4207:     PetscCall(MatCreate(PETSC_COMM_SELF, &newmat));
4208:     PetscCall(MatSetSizes(newmat, am * bm, an * bn, am * bm, an * bn));
4209:     PetscCall(MatSetType(newmat, MATAIJ));
4210:     i[0] = 0;
4211:     for (m = 0; m < am; ++m) {
4212:       for (p = 0; p < bm; ++p) {
4213:         i[m * bm + p + 1] = i[m * bm + p] + (a->i[m + 1] - a->i[m]) * (b->i[p + 1] - b->i[p]);
4214:         for (n = a->i[m]; n < a->i[m + 1]; ++n) {
4215:           for (q = b->i[p]; q < b->i[p + 1]; ++q) j[nnz++] = a->j[n] * bn + b->j[q];
4216:         }
4217:       }
4218:     }
4219:     PetscCall(MatSeqAIJSetPreallocationCSR(newmat, i, j, NULL));
4220:     *C = newmat;
4221:     PetscCall(PetscFree2(i, j));
4222:     nnz = 0;
4223:   }
4224:   PetscCall(MatSeqAIJGetArray(*C, &v));
4225:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
4226:   PetscCall(MatSeqAIJGetArrayRead(B, &ba));
4227:   for (m = 0; m < am; ++m) {
4228:     for (p = 0; p < bm; ++p) {
4229:       for (n = a->i[m]; n < a->i[m + 1]; ++n) {
4230:         for (q = b->i[p]; q < b->i[p + 1]; ++q) v[nnz++] = aa[n] * ba[q];
4231:       }
4232:     }
4233:   }
4234:   PetscCall(MatSeqAIJRestoreArray(*C, &v));
4235:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
4236:   PetscCall(MatSeqAIJRestoreArrayRead(B, &ba));
4237:   PetscFunctionReturn(PETSC_SUCCESS);
4238: }

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

4243: /*
4244:     Computes (B'*A')' since computing B*A directly is untenable

4246:                n                       p                          p
4247:         [             ]       [             ]         [                 ]
4248:       m [      A      ]  *  n [       B     ]   =   m [         C       ]
4249:         [             ]       [             ]         [                 ]

4251: */
4252: PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat A, Mat B, Mat C)
4253: {
4254:   Mat_SeqDense      *sub_a = (Mat_SeqDense *)A->data;
4255:   Mat_SeqAIJ        *sub_b = (Mat_SeqAIJ *)B->data;
4256:   Mat_SeqDense      *sub_c = (Mat_SeqDense *)C->data;
4257:   PetscInt           i, j, n, m, q, p;
4258:   const PetscInt    *ii, *idx;
4259:   const PetscScalar *b, *a, *a_q;
4260:   PetscScalar       *c, *c_q;
4261:   PetscInt           clda = sub_c->lda;
4262:   PetscInt           alda = sub_a->lda;

4264:   PetscFunctionBegin;
4265:   m = A->rmap->n;
4266:   n = A->cmap->n;
4267:   p = B->cmap->n;
4268:   a = sub_a->v;
4269:   b = sub_b->a;
4270:   c = sub_c->v;
4271:   if (clda == m) {
4272:     PetscCall(PetscArrayzero(c, m * p));
4273:   } else {
4274:     for (j = 0; j < p; j++)
4275:       for (i = 0; i < m; i++) c[j * clda + i] = 0.0;
4276:   }
4277:   ii  = sub_b->i;
4278:   idx = sub_b->j;
4279:   for (i = 0; i < n; i++) {
4280:     q = ii[i + 1] - ii[i];
4281:     while (q-- > 0) {
4282:       c_q = c + clda * (*idx);
4283:       a_q = a + alda * i;
4284:       PetscKernelAXPY(c_q, *b, a_q, m);
4285:       idx++;
4286:       b++;
4287:     }
4288:   }
4289:   PetscFunctionReturn(PETSC_SUCCESS);
4290: }

4292: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat A, Mat B, PetscReal fill, Mat C)
4293: {
4294:   PetscInt  m = A->rmap->n, n = B->cmap->n;
4295:   PetscBool cisdense;

4297:   PetscFunctionBegin;
4298:   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);
4299:   PetscCall(MatSetSizes(C, m, n, m, n));
4300:   PetscCall(MatSetBlockSizesFromMats(C, A, B));
4301:   PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATSEQDENSECUDA, MATSEQDENSEHIP, ""));
4302:   if (!cisdense) PetscCall(MatSetType(C, MATDENSE));
4303:   PetscCall(MatSetUp(C));

4305:   C->ops->matmultnumeric = MatMatMultNumeric_SeqDense_SeqAIJ;
4306:   PetscFunctionReturn(PETSC_SUCCESS);
4307: }

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

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

4316:    Level: beginner

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

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

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

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

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

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

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

4344:   Level: beginner

4346:    Note:
4347:    Subclasses include `MATAIJCUSPARSE`, `MATAIJPERM`, `MATAIJSELL`, `MATAIJMKL`, `MATAIJCRL`, and also automatically switches over to use inodes when
4348:    enough exist.

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

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

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

4359:   Level: beginner

4361:    Note:
4362:    This matrix type is identical to `MATSEQAIJCRL` when constructed with a single process communicator,
4363:    and `MATMPIAIJCRL` otherwise.  As a result, for single process communicators,
4364:    `MatSeqAIJSetPreallocation()` is supported, and similarly `MatMPIAIJSetPreallocation()` is supported
4365:    for communicators controlling multiple processes.  It is recommended that you call both of
4366:    the above preallocation routines for simplicity.

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

4371: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *);
4372: #if defined(PETSC_HAVE_ELEMENTAL)
4373: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
4374: #endif
4375: #if defined(PETSC_HAVE_SCALAPACK)
4376: PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
4377: #endif
4378: #if defined(PETSC_HAVE_HYPRE)
4379: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType, MatReuse, Mat *);
4380: #endif

4382: PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat, MatType, MatReuse, Mat *);
4383: PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat, MatType, MatReuse, Mat *);
4384: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_IS_XAIJ(Mat);

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

4389:   Not Collective

4391:   Input Parameter:
4392: . A - a `MATSEQAIJ` matrix

4394:   Output Parameter:
4395: . array - pointer to the data

4397:   Level: intermediate

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

4402: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRestoreArray()`, `MatSeqAIJGetArrayF90()`
4403: @*/
4404: PetscErrorCode MatSeqAIJGetArray(Mat A, PetscScalar *array[])
4405: {
4406:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4408:   PetscFunctionBegin;
4409:   if (aij->ops->getarray) {
4410:     PetscCall((*aij->ops->getarray)(A, array));
4411:   } else {
4412:     *array = aij->a;
4413:   }
4414:   PetscFunctionReturn(PETSC_SUCCESS);
4415: }

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

4420:   Not Collective

4422:   Input Parameters:
4423: + A     - a `MATSEQAIJ` matrix
4424: - array - pointer to the data

4426:   Level: intermediate

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

4431: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayF90()`
4432: @*/
4433: PetscErrorCode MatSeqAIJRestoreArray(Mat A, PetscScalar *array[])
4434: {
4435:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4437:   PetscFunctionBegin;
4438:   if (aij->ops->restorearray) {
4439:     PetscCall((*aij->ops->restorearray)(A, array));
4440:   } else {
4441:     *array = NULL;
4442:   }
4443:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
4444:   PetscCall(PetscObjectStateIncrease((PetscObject)A));
4445:   PetscFunctionReturn(PETSC_SUCCESS);
4446: }

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

4451:   Not Collective; No Fortran Support

4453:   Input Parameter:
4454: . A - a `MATSEQAIJ` matrix

4456:   Output Parameter:
4457: . array - pointer to the data

4459:   Level: intermediate

4461: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayRead()`
4462: @*/
4463: PetscErrorCode MatSeqAIJGetArrayRead(Mat A, const PetscScalar *array[])
4464: {
4465:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4467:   PetscFunctionBegin;
4468:   if (aij->ops->getarrayread) {
4469:     PetscCall((*aij->ops->getarrayread)(A, array));
4470:   } else {
4471:     *array = aij->a;
4472:   }
4473:   PetscFunctionReturn(PETSC_SUCCESS);
4474: }

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

4479:   Not Collective; No Fortran Support

4481:   Input Parameter:
4482: . A - a `MATSEQAIJ` matrix

4484:   Output Parameter:
4485: . array - pointer to the data

4487:   Level: intermediate

4489: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()`
4490: @*/
4491: PetscErrorCode MatSeqAIJRestoreArrayRead(Mat A, const PetscScalar *array[])
4492: {
4493:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4495:   PetscFunctionBegin;
4496:   if (aij->ops->restorearrayread) {
4497:     PetscCall((*aij->ops->restorearrayread)(A, array));
4498:   } else {
4499:     *array = NULL;
4500:   }
4501:   PetscFunctionReturn(PETSC_SUCCESS);
4502: }

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

4507:   Not Collective; No Fortran Support

4509:   Input Parameter:
4510: . A - a `MATSEQAIJ` matrix

4512:   Output Parameter:
4513: . array - pointer to the data

4515:   Level: intermediate

4517: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArrayRead()`
4518: @*/
4519: PetscErrorCode MatSeqAIJGetArrayWrite(Mat A, PetscScalar *array[])
4520: {
4521:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4523:   PetscFunctionBegin;
4524:   if (aij->ops->getarraywrite) {
4525:     PetscCall((*aij->ops->getarraywrite)(A, array));
4526:   } else {
4527:     *array = aij->a;
4528:   }
4529:   PetscCall(MatSeqAIJInvalidateDiagonal(A));
4530:   PetscCall(PetscObjectStateIncrease((PetscObject)A));
4531:   PetscFunctionReturn(PETSC_SUCCESS);
4532: }

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

4537:   Not Collective; No Fortran Support

4539:   Input Parameter:
4540: . A - a MATSEQAIJ matrix

4542:   Output Parameter:
4543: . array - pointer to the data

4545:   Level: intermediate

4547: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()`
4548: @*/
4549: PetscErrorCode MatSeqAIJRestoreArrayWrite(Mat A, PetscScalar *array[])
4550: {
4551:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4553:   PetscFunctionBegin;
4554:   if (aij->ops->restorearraywrite) {
4555:     PetscCall((*aij->ops->restorearraywrite)(A, array));
4556:   } else {
4557:     *array = NULL;
4558:   }
4559:   PetscFunctionReturn(PETSC_SUCCESS);
4560: }

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

4565:   Not Collective; No Fortran Support

4567:   Input Parameter:
4568: . mat - a matrix of type `MATSEQAIJ` or its subclasses

4570:   Output Parameters:
4571: + i     - row map array of the matrix
4572: . j     - column index array of the matrix
4573: . a     - data array of the matrix
4574: - mtype - memory type of the arrays

4576:   Level: developer

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

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

4585: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArray()`, `MatSeqAIJGetArrayRead()`
4586: @*/
4587: PetscErrorCode MatSeqAIJGetCSRAndMemType(Mat mat, const PetscInt *i[], const PetscInt *j[], PetscScalar *a[], PetscMemType *mtype)
4588: {
4589:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;

4591:   PetscFunctionBegin;
4592:   PetscCheck(mat->preallocated, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "matrix is not preallocated");
4593:   if (aij->ops->getcsrandmemtype) {
4594:     PetscCall((*aij->ops->getcsrandmemtype)(mat, i, j, a, mtype));
4595:   } else {
4596:     if (i) *i = aij->i;
4597:     if (j) *j = aij->j;
4598:     if (a) *a = aij->a;
4599:     if (mtype) *mtype = PETSC_MEMTYPE_HOST;
4600:   }
4601:   PetscFunctionReturn(PETSC_SUCCESS);
4602: }

4604: /*@
4605:   MatSeqAIJGetMaxRowNonzeros - returns the maximum number of nonzeros in any row

4607:   Not Collective

4609:   Input Parameter:
4610: . A - a `MATSEQAIJ` matrix

4612:   Output Parameter:
4613: . nz - the maximum number of nonzeros in any row

4615:   Level: intermediate

4617: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRestoreArray()`, `MatSeqAIJGetArrayF90()`
4618: @*/
4619: PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat A, PetscInt *nz)
4620: {
4621:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;

4623:   PetscFunctionBegin;
4624:   *nz = aij->rmax;
4625:   PetscFunctionReturn(PETSC_SUCCESS);
4626: }

4628: static PetscErrorCode MatCOOStructDestroy_SeqAIJ(void *data)
4629: {
4630:   MatCOOStruct_SeqAIJ *coo = (MatCOOStruct_SeqAIJ *)data;

4632:   PetscFunctionBegin;
4633:   PetscCall(PetscFree(coo->perm));
4634:   PetscCall(PetscFree(coo->jmap));
4635:   PetscCall(PetscFree(coo));
4636:   PetscFunctionReturn(PETSC_SUCCESS);
4637: }

4639: PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat mat, PetscCount coo_n, PetscInt coo_i[], PetscInt coo_j[])
4640: {
4641:   MPI_Comm             comm;
4642:   PetscInt            *i, *j;
4643:   PetscInt             M, N, row, iprev;
4644:   PetscCount           k, p, q, nneg, nnz, start, end; /* Index the coo array, so use PetscCount as their type */
4645:   PetscInt            *Ai;                             /* Change to PetscCount once we use it for row pointers */
4646:   PetscInt            *Aj;
4647:   PetscScalar         *Aa;
4648:   Mat_SeqAIJ          *seqaij = (Mat_SeqAIJ *)mat->data;
4649:   MatType              rtype;
4650:   PetscCount          *perm, *jmap;
4651:   MatCOOStruct_SeqAIJ *coo;
4652:   PetscBool            isorted;
4653:   PetscBool            hypre;
4654:   const char          *name;

4656:   PetscFunctionBegin;
4657:   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
4658:   PetscCall(MatGetSize(mat, &M, &N));
4659:   i = coo_i;
4660:   j = coo_j;
4661:   PetscCall(PetscMalloc1(coo_n, &perm));

4663:   /* 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) */
4664:   isorted = PETSC_TRUE;
4665:   iprev   = PETSC_INT_MIN;
4666:   for (k = 0; k < coo_n; k++) {
4667:     if (j[k] < 0) i[k] = -1;
4668:     if (isorted) {
4669:       if (i[k] < iprev) isorted = PETSC_FALSE;
4670:       else iprev = i[k];
4671:     }
4672:     perm[k] = k;
4673:   }

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

4678:   /* Advance k to the first row with a non-negative index */
4679:   for (k = 0; k < coo_n; k++)
4680:     if (i[k] >= 0) break;
4681:   nneg = k;
4682:   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 */
4683:   nnz = 0;                                          /* Total number of unique nonzeros to be counted */
4684:   jmap++;                                           /* Inc jmap by 1 for convenience */

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

4690:   PetscCall(PetscObjectGetName((PetscObject)mat, &name));
4691:   PetscCall(PetscStrcmp("_internal_COO_mat_for_hypre", name, &hypre));

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

4700:     /* get [start,end) indices for this row; also check if cols in this row are strictly sorted */
4701:     row             = i[k];
4702:     start           = k;
4703:     jprev           = PETSC_INT_MIN;
4704:     strictly_sorted = PETSC_TRUE;
4705:     while (k < coo_n && i[k] == row) {
4706:       if (strictly_sorted) {
4707:         if (j[k] <= jprev) strictly_sorted = PETSC_FALSE;
4708:         else jprev = j[k];
4709:       }
4710:       k++;
4711:     }
4712:     end = k;

4714:     /* hack for HYPRE: swap min column to diag so that diagonal values will go first */
4715:     if (hypre) {
4716:       PetscInt  minj    = PETSC_INT_MAX;
4717:       PetscBool hasdiag = PETSC_FALSE;

4719:       if (strictly_sorted) { // fast path to swap the first and the diag
4720:         PetscCount tmp;
4721:         for (p = start; p < end; p++) {
4722:           if (j[p] == row && p != start) {
4723:             j[p]        = j[start];
4724:             j[start]    = row;
4725:             tmp         = perm[start];
4726:             perm[start] = perm[p];
4727:             perm[p]     = tmp;
4728:             break;
4729:           }
4730:         }
4731:       } else {
4732:         for (p = start; p < end; p++) {
4733:           hasdiag = (PetscBool)(hasdiag || (j[p] == row));
4734:           minj    = PetscMin(minj, j[p]);
4735:         }

4737:         if (hasdiag) {
4738:           for (p = start; p < end; p++) {
4739:             if (j[p] == minj) j[p] = row;
4740:             else if (j[p] == row) j[p] = minj;
4741:           }
4742:         }
4743:       }
4744:     }
4745:     // sort by columns in a row
4746:     if (!strictly_sorted) PetscCall(PetscSortIntWithCountArray(end - start, j + start, perm + start));

4748:     if (strictly_sorted) { // fast path to set Aj[], jmap[], Ai[], nnz, q
4749:       for (p = start; p < end; p++, q++) {
4750:         Aj[q]   = j[p];
4751:         jmap[q] = 1;
4752:       }
4753:       PetscCall(PetscIntCast(end - start, Ai + row));
4754:       nnz += Ai[row]; // q is already advanced
4755:     } else {
4756:       /* Find number of unique col entries in this row */
4757:       Aj[q]   = j[start]; /* Log the first nonzero in this row */
4758:       jmap[q] = 1;        /* Number of repeats of this nonzero entry */
4759:       Ai[row] = 1;
4760:       nnz++;

4762:       for (p = start + 1; p < end; p++) { /* Scan remaining nonzero in this row */
4763:         if (j[p] != j[p - 1]) {           /* Meet a new nonzero */
4764:           q++;
4765:           jmap[q] = 1;
4766:           Aj[q]   = j[p];
4767:           Ai[row]++;
4768:           nnz++;
4769:         } else {
4770:           jmap[q]++;
4771:         }
4772:       }
4773:       q++; /* Move to next row and thus next unique nonzero */
4774:     }
4775:   }

4777:   Ai--; /* Back to the beginning of Ai[] */
4778:   for (k = 0; k < M; k++) Ai[k + 1] += Ai[k];
4779:   jmap--; // Back to the beginning of jmap[]
4780:   jmap[0] = 0;
4781:   for (k = 0; k < nnz; k++) jmap[k + 1] += jmap[k];

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

4787:     PetscCall(PetscMalloc1(nnz + 1, &jmap_new));
4788:     PetscCall(PetscArraycpy(jmap_new, jmap, nnz + 1));
4789:     PetscCall(PetscFree(jmap));
4790:     jmap = jmap_new;

4792:     PetscCall(PetscShmgetAllocateArray(nnz, sizeof(PetscInt), (void **)&Aj_new));
4793:     PetscCall(PetscArraycpy(Aj_new, Aj, nnz));
4794:     PetscCall(PetscShmgetDeallocateArray((void **)&Aj));
4795:     Aj = Aj_new;
4796:   }

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

4801:     PetscCall(PetscMalloc1(coo_n - nneg, &perm_new));
4802:     PetscCall(PetscArraycpy(perm_new, perm + nneg, coo_n - nneg));
4803:     PetscCall(PetscFree(perm));
4804:     perm = perm_new;
4805:   }

4807:   PetscCall(MatGetRootType_Private(mat, &rtype));
4808:   PetscCall(PetscShmgetAllocateArray(nnz, sizeof(PetscScalar), (void **)&Aa));
4809:   PetscCall(PetscArrayzero(Aa, nnz));
4810:   PetscCall(MatSetSeqAIJWithArrays_private(PETSC_COMM_SELF, M, N, Ai, Aj, Aa, rtype, mat));

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

4814:   // Put the COO struct in a container and then attach that to the matrix
4815:   PetscCall(PetscMalloc1(1, &coo));
4816:   PetscCall(PetscIntCast(nnz, &coo->nz));
4817:   coo->n    = coo_n;
4818:   coo->Atot = coo_n - nneg; // Annz is seqaij->nz, so no need to record that again
4819:   coo->jmap = jmap;         // of length nnz+1
4820:   coo->perm = perm;
4821:   PetscCall(PetscObjectContainerCompose((PetscObject)mat, "__PETSc_MatCOOStruct_Host", coo, MatCOOStructDestroy_SeqAIJ));
4822:   PetscFunctionReturn(PETSC_SUCCESS);
4823: }

4825: static PetscErrorCode MatSetValuesCOO_SeqAIJ(Mat A, const PetscScalar v[], InsertMode imode)
4826: {
4827:   Mat_SeqAIJ          *aseq = (Mat_SeqAIJ *)A->data;
4828:   PetscCount           i, j, Annz = aseq->nz;
4829:   PetscCount          *perm, *jmap;
4830:   PetscScalar         *Aa;
4831:   PetscContainer       container;
4832:   MatCOOStruct_SeqAIJ *coo;

4834:   PetscFunctionBegin;
4835:   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_MatCOOStruct_Host", (PetscObject *)&container));
4836:   PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not found MatCOOStruct on this matrix");
4837:   PetscCall(PetscContainerGetPointer(container, (void **)&coo));
4838:   perm = coo->perm;
4839:   jmap = coo->jmap;
4840:   PetscCall(MatSeqAIJGetArray(A, &Aa));
4841:   for (i = 0; i < Annz; i++) {
4842:     PetscScalar sum = 0.0;
4843:     for (j = jmap[i]; j < jmap[i + 1]; j++) sum += v[perm[j]];
4844:     Aa[i] = (imode == INSERT_VALUES ? 0.0 : Aa[i]) + sum;
4845:   }
4846:   PetscCall(MatSeqAIJRestoreArray(A, &Aa));
4847:   PetscFunctionReturn(PETSC_SUCCESS);
4848: }

4850: #if defined(PETSC_HAVE_CUDA)
4851: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCUSPARSE(Mat, MatType, MatReuse, Mat *);
4852: #endif
4853: #if defined(PETSC_HAVE_HIP)
4854: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJHIPSPARSE(Mat, MatType, MatReuse, Mat *);
4855: #endif
4856: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
4857: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJKokkos(Mat, MatType, MatReuse, Mat *);
4858: #endif

4860: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat B)
4861: {
4862:   Mat_SeqAIJ *b;
4863:   PetscMPIInt size;

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

4869:   PetscCall(PetscNew(&b));

4871:   B->data   = (void *)b;
4872:   B->ops[0] = MatOps_Values;
4873:   if (B->sortedfull) B->ops->setvalues = MatSetValues_SeqAIJ_SortedFull;

4875:   b->row                = NULL;
4876:   b->col                = NULL;
4877:   b->icol               = NULL;
4878:   b->reallocs           = 0;
4879:   b->ignorezeroentries  = PETSC_FALSE;
4880:   b->roworiented        = PETSC_TRUE;
4881:   b->nonew              = 0;
4882:   b->diag               = NULL;
4883:   b->solve_work         = NULL;
4884:   B->spptr              = NULL;
4885:   b->saved_values       = NULL;
4886:   b->idiag              = NULL;
4887:   b->mdiag              = NULL;
4888:   b->ssor_work          = NULL;
4889:   b->omega              = 1.0;
4890:   b->fshift             = 0.0;
4891:   b->idiagvalid         = PETSC_FALSE;
4892:   b->ibdiagvalid        = PETSC_FALSE;
4893:   b->keepnonzeropattern = PETSC_FALSE;

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

4955: /*
4956:     Given a matrix generated with MatGetFactor() duplicates all the information in A into C
4957: */
4958: PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat C, Mat A, MatDuplicateOption cpvalues, PetscBool mallocmatspace)
4959: {
4960:   Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data, *a = (Mat_SeqAIJ *)A->data;
4961:   PetscInt    m = A->rmap->n, i;

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

4966:   C->factortype    = A->factortype;
4967:   c->row           = NULL;
4968:   c->col           = NULL;
4969:   c->icol          = NULL;
4970:   c->reallocs      = 0;
4971:   c->diagonaldense = a->diagonaldense;

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

4975:   if (A->preallocated) {
4976:     PetscCall(PetscLayoutReference(A->rmap, &C->rmap));
4977:     PetscCall(PetscLayoutReference(A->cmap, &C->cmap));

4979:     if (!A->hash_active) {
4980:       PetscCall(PetscMalloc1(m, &c->imax));
4981:       PetscCall(PetscMemcpy(c->imax, a->imax, m * sizeof(PetscInt)));
4982:       PetscCall(PetscMalloc1(m, &c->ilen));
4983:       PetscCall(PetscMemcpy(c->ilen, a->ilen, m * sizeof(PetscInt)));

4985:       /* allocate the matrix space */
4986:       if (mallocmatspace) {
4987:         PetscCall(PetscShmgetAllocateArray(a->i[m], sizeof(PetscScalar), (void **)&c->a));
4988:         PetscCall(PetscShmgetAllocateArray(a->i[m], sizeof(PetscInt), (void **)&c->j));
4989:         PetscCall(PetscShmgetAllocateArray(m + 1, sizeof(PetscInt), (void **)&c->i));
4990:         PetscCall(PetscArraycpy(c->i, a->i, m + 1));
4991:         c->free_a  = PETSC_TRUE;
4992:         c->free_ij = PETSC_TRUE;
4993:         if (m > 0) {
4994:           PetscCall(PetscArraycpy(c->j, a->j, a->i[m]));
4995:           if (cpvalues == MAT_COPY_VALUES) {
4996:             const PetscScalar *aa;

4998:             PetscCall(MatSeqAIJGetArrayRead(A, &aa));
4999:             PetscCall(PetscArraycpy(c->a, aa, a->i[m]));
5000:             PetscCall(MatSeqAIJGetArrayRead(A, &aa));
5001:           } else {
5002:             PetscCall(PetscArrayzero(c->a, a->i[m]));
5003:           }
5004:         }
5005:       }
5006:       C->preallocated = PETSC_TRUE;
5007:     } else {
5008:       PetscCheck(mallocmatspace, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Cannot malloc matrix memory from a non-preallocated matrix");
5009:       PetscCall(MatSetUp(C));
5010:     }

5012:     c->ignorezeroentries = a->ignorezeroentries;
5013:     c->roworiented       = a->roworiented;
5014:     c->nonew             = a->nonew;
5015:     if (a->diag) {
5016:       PetscCall(PetscMalloc1(m + 1, &c->diag));
5017:       PetscCall(PetscMemcpy(c->diag, a->diag, m * sizeof(PetscInt)));
5018:     } else c->diag = NULL;

5020:     c->solve_work         = NULL;
5021:     c->saved_values       = NULL;
5022:     c->idiag              = NULL;
5023:     c->ssor_work          = NULL;
5024:     c->keepnonzeropattern = a->keepnonzeropattern;

5026:     c->rmax  = a->rmax;
5027:     c->nz    = a->nz;
5028:     c->maxnz = a->nz; /* Since we allocate exactly the right amount */

5030:     c->compressedrow.use   = a->compressedrow.use;
5031:     c->compressedrow.nrows = a->compressedrow.nrows;
5032:     if (a->compressedrow.use) {
5033:       i = a->compressedrow.nrows;
5034:       PetscCall(PetscMalloc2(i + 1, &c->compressedrow.i, i, &c->compressedrow.rindex));
5035:       PetscCall(PetscArraycpy(c->compressedrow.i, a->compressedrow.i, i + 1));
5036:       PetscCall(PetscArraycpy(c->compressedrow.rindex, a->compressedrow.rindex, i));
5037:     } else {
5038:       c->compressedrow.use    = PETSC_FALSE;
5039:       c->compressedrow.i      = NULL;
5040:       c->compressedrow.rindex = NULL;
5041:     }
5042:     c->nonzerorowcnt = a->nonzerorowcnt;
5043:     C->nonzerostate  = A->nonzerostate;

5045:     PetscCall(MatDuplicate_SeqAIJ_Inode(A, cpvalues, &C));
5046:   }
5047:   PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist, &((PetscObject)C)->qlist));
5048:   PetscFunctionReturn(PETSC_SUCCESS);
5049: }

5051: PetscErrorCode MatDuplicate_SeqAIJ(Mat A, MatDuplicateOption cpvalues, Mat *B)
5052: {
5053:   PetscFunctionBegin;
5054:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
5055:   PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n));
5056:   if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) PetscCall(MatSetBlockSizesFromMats(*B, A, A));
5057:   PetscCall(MatSetType(*B, ((PetscObject)A)->type_name));
5058:   PetscCall(MatDuplicateNoCreate_SeqAIJ(*B, A, cpvalues, PETSC_TRUE));
5059:   PetscFunctionReturn(PETSC_SUCCESS);
5060: }

5062: PetscErrorCode MatLoad_SeqAIJ(Mat newMat, PetscViewer viewer)
5063: {
5064:   PetscBool isbinary, ishdf5;

5066:   PetscFunctionBegin;
5069:   /* force binary viewer to load .info file if it has not yet done so */
5070:   PetscCall(PetscViewerSetUp(viewer));
5071:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
5072:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
5073:   if (isbinary) {
5074:     PetscCall(MatLoad_SeqAIJ_Binary(newMat, viewer));
5075:   } else if (ishdf5) {
5076: #if defined(PETSC_HAVE_HDF5)
5077:     PetscCall(MatLoad_AIJ_HDF5(newMat, viewer));
5078: #else
5079:     SETERRQ(PetscObjectComm((PetscObject)newMat), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
5080: #endif
5081:   } else {
5082:     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);
5083:   }
5084:   PetscFunctionReturn(PETSC_SUCCESS);
5085: }

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

5092:   PetscFunctionBegin;
5093:   PetscCall(PetscViewerSetUp(viewer));

5095:   /* read in matrix header */
5096:   PetscCall(PetscViewerBinaryRead(viewer, header, 4, NULL, PETSC_INT));
5097:   PetscCheck(header[0] == MAT_FILE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a matrix object in file");
5098:   M  = header[1];
5099:   N  = header[2];
5100:   nz = header[3];
5101:   PetscCheck(M >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Matrix row size (%" PetscInt_FMT ") in file is negative", M);
5102:   PetscCheck(N >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Matrix column size (%" PetscInt_FMT ") in file is negative", N);
5103:   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Matrix stored in special format on disk, cannot load as SeqAIJ");

5105:   /* set block sizes from the viewer's .info file */
5106:   PetscCall(MatLoad_Binary_BlockSizes(mat, viewer));
5107:   /* set local and global sizes if not set already */
5108:   if (mat->rmap->n < 0) mat->rmap->n = M;
5109:   if (mat->cmap->n < 0) mat->cmap->n = N;
5110:   if (mat->rmap->N < 0) mat->rmap->N = M;
5111:   if (mat->cmap->N < 0) mat->cmap->N = N;
5112:   PetscCall(PetscLayoutSetUp(mat->rmap));
5113:   PetscCall(PetscLayoutSetUp(mat->cmap));

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

5119:   /* read in row lengths */
5120:   PetscCall(PetscMalloc1(M, &rowlens));
5121:   PetscCall(PetscViewerBinaryRead(viewer, rowlens, M, NULL, PETSC_INT));
5122:   /* check if sum(rowlens) is same as nz */
5123:   sum = 0;
5124:   for (i = 0; i < M; i++) sum += rowlens[i];
5125:   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);
5126:   /* preallocate and check sizes */
5127:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(mat, 0, rowlens));
5128:   PetscCall(MatGetSize(mat, &rows, &cols));
5129:   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);
5130:   /* store row lengths */
5131:   PetscCall(PetscArraycpy(a->ilen, rowlens, M));
5132:   PetscCall(PetscFree(rowlens));

5134:   /* fill in "i" row pointers */
5135:   a->i[0] = 0;
5136:   for (i = 0; i < M; i++) a->i[i + 1] = a->i[i] + a->ilen[i];
5137:   /* read in "j" column indices */
5138:   PetscCall(PetscViewerBinaryRead(viewer, a->j, nz, NULL, PETSC_INT));
5139:   /* read in "a" nonzero values */
5140:   PetscCall(PetscViewerBinaryRead(viewer, a->a, nz, NULL, PETSC_SCALAR));

5142:   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
5143:   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
5144:   PetscFunctionReturn(PETSC_SUCCESS);
5145: }

5147: PetscErrorCode MatEqual_SeqAIJ(Mat A, Mat B, PetscBool *flg)
5148: {
5149:   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
5150:   const PetscScalar *aa, *ba;
5151: #if defined(PETSC_USE_COMPLEX)
5152:   PetscInt k;
5153: #endif

5155:   PetscFunctionBegin;
5156:   /* If the  matrix dimensions are not equal,or no of nonzeros */
5157:   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) || (a->nz != b->nz)) {
5158:     *flg = PETSC_FALSE;
5159:     PetscFunctionReturn(PETSC_SUCCESS);
5160:   }

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

5166:   /* if a->j are the same */
5167:   PetscCall(PetscArraycmp(a->j, b->j, a->nz, flg));
5168:   if (!*flg) PetscFunctionReturn(PETSC_SUCCESS);

5170:   PetscCall(MatSeqAIJGetArrayRead(A, &aa));
5171:   PetscCall(MatSeqAIJGetArrayRead(B, &ba));
5172:   /* if a->a are the same */
5173: #if defined(PETSC_USE_COMPLEX)
5174:   for (k = 0; k < a->nz; k++) {
5175:     if (PetscRealPart(aa[k]) != PetscRealPart(ba[k]) || PetscImaginaryPart(aa[k]) != PetscImaginaryPart(ba[k])) {
5176:       *flg = PETSC_FALSE;
5177:       PetscFunctionReturn(PETSC_SUCCESS);
5178:     }
5179:   }
5180: #else
5181:   PetscCall(PetscArraycmp(aa, ba, a->nz, flg));
5182: #endif
5183:   PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
5184:   PetscCall(MatSeqAIJRestoreArrayRead(B, &ba));
5185:   PetscFunctionReturn(PETSC_SUCCESS);
5186: }

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

5192:   Collective

5194:   Input Parameters:
5195: + comm - must be an MPI communicator of size 1
5196: . m    - number of rows
5197: . n    - number of columns
5198: . i    - row indices; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
5199: . j    - column indices
5200: - a    - matrix values

5202:   Output Parameter:
5203: . mat - the matrix

5205:   Level: intermediate

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

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

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

5215:   The format which is used for the sparse matrix input, is equivalent to a
5216:   row-major ordering.. i.e for the following matrix, the input data expected is
5217:   as shown
5218: .vb
5219:         1 0 0
5220:         2 0 3
5221:         4 5 6

5223:         i =  {0,1,3,6}  [size = nrow+1  = 3+1]
5224:         j =  {0,0,2,0,1,2}  [size = 6]; values must be sorted for each row
5225:         v =  {1,2,3,4,5,6}  [size = 6]
5226: .ve

5228: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateAIJ()`, `MatCreateSeqAIJ()`, `MatCreateMPIAIJWithArrays()`, `MatMPIAIJSetPreallocationCSR()`
5229: @*/
5230: PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt i[], PetscInt j[], PetscScalar a[], Mat *mat)
5231: {
5232:   PetscInt    ii;
5233:   Mat_SeqAIJ *aij;
5234:   PetscInt    jj;

5236:   PetscFunctionBegin;
5237:   PetscCheck(m <= 0 || i[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "i (row indices) must start with 0");
5238:   PetscCall(MatCreate(comm, mat));
5239:   PetscCall(MatSetSizes(*mat, m, n, m, n));
5240:   /* PetscCall(MatSetBlockSizes(*mat,,)); */
5241:   PetscCall(MatSetType(*mat, MATSEQAIJ));
5242:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*mat, MAT_SKIP_ALLOCATION, NULL));
5243:   aij = (Mat_SeqAIJ *)(*mat)->data;
5244:   PetscCall(PetscMalloc1(m, &aij->imax));
5245:   PetscCall(PetscMalloc1(m, &aij->ilen));

5247:   aij->i       = i;
5248:   aij->j       = j;
5249:   aij->a       = a;
5250:   aij->nonew   = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
5251:   aij->free_a  = PETSC_FALSE;
5252:   aij->free_ij = PETSC_FALSE;

5254:   for (ii = 0, aij->nonzerorowcnt = 0, aij->rmax = 0; ii < m; ii++) {
5255:     aij->ilen[ii] = aij->imax[ii] = i[ii + 1] - i[ii];
5256:     if (PetscDefined(USE_DEBUG)) {
5257:       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]);
5258:       for (jj = i[ii] + 1; jj < i[ii + 1]; jj++) {
5259:         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);
5260:         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);
5261:       }
5262:     }
5263:   }
5264:   if (PetscDefined(USE_DEBUG)) {
5265:     for (ii = 0; ii < aij->i[m]; ii++) {
5266:       PetscCheck(j[ii] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative column index at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]);
5267:       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);
5268:     }
5269:   }

5271:   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
5272:   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
5273:   PetscFunctionReturn(PETSC_SUCCESS);
5274: }

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

5280:   Collective

5282:   Input Parameters:
5283: + comm - must be an MPI communicator of size 1
5284: . m    - number of rows
5285: . n    - number of columns
5286: . i    - row indices
5287: . j    - column indices
5288: . a    - matrix values
5289: . nz   - number of nonzeros
5290: - idx  - if the `i` and `j` indices start with 1 use `PETSC_TRUE` otherwise use `PETSC_FALSE`

5292:   Output Parameter:
5293: . mat - the matrix

5295:   Level: intermediate

5297:   Example:
5298:   For the following matrix, the input data expected is as shown (using 0 based indexing)
5299: .vb
5300:         1 0 0
5301:         2 0 3
5302:         4 5 6

5304:         i =  {0,1,1,2,2,2}
5305:         j =  {0,0,2,0,1,2}
5306:         v =  {1,2,3,4,5,6}
5307: .ve

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

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

5319:   PetscFunctionBegin;
5320:   PetscCall(PetscCalloc1(m, &nnz));
5321:   for (ii = 0; ii < nz; ii++) nnz[i[ii] - !!idx] += 1;
5322:   PetscCall(MatCreate(comm, mat));
5323:   PetscCall(MatSetSizes(*mat, m, n, m, n));
5324:   PetscCall(MatSetType(*mat, MATSEQAIJ));
5325:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*mat, 0, nnz));
5326:   for (ii = 0; ii < nz; ii++) {
5327:     if (idx) {
5328:       row = i[ii] - 1;
5329:       col = j[ii] - 1;
5330:     } else {
5331:       row = i[ii];
5332:       col = j[ii];
5333:     }
5334:     PetscCall(MatSetValues(*mat, one, &row, one, &col, &a[ii], ADD_VALUES));
5335:   }
5336:   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
5337:   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
5338:   PetscCall(PetscFree(nnz));
5339:   PetscFunctionReturn(PETSC_SUCCESS);
5340: }

5342: PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat A)
5343: {
5344:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

5346:   PetscFunctionBegin;
5347:   a->idiagvalid  = PETSC_FALSE;
5348:   a->ibdiagvalid = PETSC_FALSE;

5350:   PetscCall(MatSeqAIJInvalidateDiagonal_Inode(A));
5351:   PetscFunctionReturn(PETSC_SUCCESS);
5352: }

5354: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm comm, Mat inmat, PetscInt n, MatReuse scall, Mat *outmat)
5355: {
5356:   PetscFunctionBegin;
5357:   PetscCall(MatCreateMPIMatConcatenateSeqMat_MPIAIJ(comm, inmat, n, scall, outmat));
5358:   PetscFunctionReturn(PETSC_SUCCESS);
5359: }

5361: /*
5362:  Permute A into C's *local* index space using rowemb,colemb.
5363:  The embedding are supposed to be injections and the above implies that the range of rowemb is a subset
5364:  of [0,m), colemb is in [0,n).
5365:  If pattern == DIFFERENT_NONZERO_PATTERN, C is preallocated according to A.
5366:  */
5367: PetscErrorCode MatSetSeqMat_SeqAIJ(Mat C, IS rowemb, IS colemb, MatStructure pattern, Mat B)
5368: {
5369:   /* If making this function public, change the error returned in this function away from _PLIB. */
5370:   Mat_SeqAIJ     *Baij;
5371:   PetscBool       seqaij;
5372:   PetscInt        m, n, *nz, i, j, count;
5373:   PetscScalar     v;
5374:   const PetscInt *rowindices, *colindices;

5376:   PetscFunctionBegin;
5377:   if (!B) PetscFunctionReturn(PETSC_SUCCESS);
5378:   /* Check to make sure the target matrix (and embeddings) are compatible with C and each other. */
5379:   PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATSEQAIJ, &seqaij));
5380:   PetscCheck(seqaij, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is of wrong type");
5381:   if (rowemb) {
5382:     PetscCall(ISGetLocalSize(rowemb, &m));
5383:     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);
5384:   } else {
5385:     PetscCheck(C->rmap->n == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is row-incompatible with the target matrix");
5386:   }
5387:   if (colemb) {
5388:     PetscCall(ISGetLocalSize(colemb, &n));
5389:     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);
5390:   } else {
5391:     PetscCheck(C->cmap->n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Input matrix is col-incompatible with the target matrix");
5392:   }

5394:   Baij = (Mat_SeqAIJ *)B->data;
5395:   if (pattern == DIFFERENT_NONZERO_PATTERN) {
5396:     PetscCall(PetscMalloc1(B->rmap->n, &nz));
5397:     for (i = 0; i < B->rmap->n; i++) nz[i] = Baij->i[i + 1] - Baij->i[i];
5398:     PetscCall(MatSeqAIJSetPreallocation(C, 0, nz));
5399:     PetscCall(PetscFree(nz));
5400:   }
5401:   if (pattern == SUBSET_NONZERO_PATTERN) PetscCall(MatZeroEntries(C));
5402:   count      = 0;
5403:   rowindices = NULL;
5404:   colindices = NULL;
5405:   if (rowemb) PetscCall(ISGetIndices(rowemb, &rowindices));
5406:   if (colemb) PetscCall(ISGetIndices(colemb, &colindices));
5407:   for (i = 0; i < B->rmap->n; i++) {
5408:     PetscInt row;
5409:     row = i;
5410:     if (rowindices) row = rowindices[i];
5411:     for (j = Baij->i[i]; j < Baij->i[i + 1]; j++) {
5412:       PetscInt col;
5413:       col = Baij->j[count];
5414:       if (colindices) col = colindices[col];
5415:       v = Baij->a[count];
5416:       PetscCall(MatSetValues(C, 1, &row, 1, &col, &v, INSERT_VALUES));
5417:       ++count;
5418:     }
5419:   }
5420:   /* FIXME: set C's nonzerostate correctly. */
5421:   /* Assembly for C is necessary. */
5422:   C->preallocated  = PETSC_TRUE;
5423:   C->assembled     = PETSC_TRUE;
5424:   C->was_assembled = PETSC_FALSE;
5425:   PetscFunctionReturn(PETSC_SUCCESS);
5426: }

5428: PetscErrorCode MatEliminateZeros_SeqAIJ(Mat A, PetscBool keep)
5429: {
5430:   Mat_SeqAIJ *a  = (Mat_SeqAIJ *)A->data;
5431:   MatScalar  *aa = a->a;
5432:   PetscInt    m = A->rmap->n, fshift = 0, fshift_prev = 0, i, k;
5433:   PetscInt   *ailen = a->ilen, *imax = a->imax, *ai = a->i, *aj = a->j, rmax = 0;

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

5471: PetscFunctionList MatSeqAIJList = NULL;

5473: /*@
5474:   MatSeqAIJSetType - Converts a `MATSEQAIJ` matrix to a subtype

5476:   Collective

5478:   Input Parameters:
5479: + mat    - the matrix object
5480: - matype - matrix type

5482:   Options Database Key:
5483: . -mat_seqaij_type  <method> - for example seqaijcrl

5485:   Level: intermediate

5487: .seealso: [](ch_matrices), `Mat`, `PCSetType()`, `VecSetType()`, `MatCreate()`, `MatType`
5488: @*/
5489: PetscErrorCode MatSeqAIJSetType(Mat mat, MatType matype)
5490: {
5491:   PetscBool sametype;
5492:   PetscErrorCode (*r)(Mat, MatType, MatReuse, Mat *);

5494:   PetscFunctionBegin;
5496:   PetscCall(PetscObjectTypeCompare((PetscObject)mat, matype, &sametype));
5497:   if (sametype) PetscFunctionReturn(PETSC_SUCCESS);

5499:   PetscCall(PetscFunctionListFind(MatSeqAIJList, matype, &r));
5500:   PetscCheck(r, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Mat type given: %s", matype);
5501:   PetscCall((*r)(mat, matype, MAT_INPLACE_MATRIX, &mat));
5502:   PetscFunctionReturn(PETSC_SUCCESS);
5503: }

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

5508:   Not Collective, No Fortran Support

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

5514:   Level: advanced

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

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

5522: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRegisterAll()`
5523: @*/
5524: PetscErrorCode MatSeqAIJRegister(const char sname[], PetscErrorCode (*function)(Mat, MatType, MatReuse, Mat *))
5525: {
5526:   PetscFunctionBegin;
5527:   PetscCall(MatInitializePackage());
5528:   PetscCall(PetscFunctionListAdd(&MatSeqAIJList, sname, function));
5529:   PetscFunctionReturn(PETSC_SUCCESS);
5530: }

5532: PetscBool MatSeqAIJRegisterAllCalled = PETSC_FALSE;

5534: /*@C
5535:   MatSeqAIJRegisterAll - Registers all of the matrix subtypes of `MATSSEQAIJ`

5537:   Not Collective

5539:   Level: advanced

5541:   Note:
5542:   This registers the versions of `MATSEQAIJ` for GPUs

5544: .seealso: [](ch_matrices), `Mat`, `MatRegisterAll()`, `MatSeqAIJRegister()`
5545: @*/
5546: PetscErrorCode MatSeqAIJRegisterAll(void)
5547: {
5548:   PetscFunctionBegin;
5549:   if (MatSeqAIJRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
5550:   MatSeqAIJRegisterAllCalled = PETSC_TRUE;

5552:   PetscCall(MatSeqAIJRegister(MATSEQAIJCRL, MatConvert_SeqAIJ_SeqAIJCRL));
5553:   PetscCall(MatSeqAIJRegister(MATSEQAIJPERM, MatConvert_SeqAIJ_SeqAIJPERM));
5554:   PetscCall(MatSeqAIJRegister(MATSEQAIJSELL, MatConvert_SeqAIJ_SeqAIJSELL));
5555: #if defined(PETSC_HAVE_MKL_SPARSE)
5556:   PetscCall(MatSeqAIJRegister(MATSEQAIJMKL, MatConvert_SeqAIJ_SeqAIJMKL));
5557: #endif
5558: #if defined(PETSC_HAVE_CUDA)
5559:   PetscCall(MatSeqAIJRegister(MATSEQAIJCUSPARSE, MatConvert_SeqAIJ_SeqAIJCUSPARSE));
5560: #endif
5561: #if defined(PETSC_HAVE_HIP)
5562:   PetscCall(MatSeqAIJRegister(MATSEQAIJHIPSPARSE, MatConvert_SeqAIJ_SeqAIJHIPSPARSE));
5563: #endif
5564: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
5565:   PetscCall(MatSeqAIJRegister(MATSEQAIJKOKKOS, MatConvert_SeqAIJ_SeqAIJKokkos));
5566: #endif
5567: #if defined(PETSC_HAVE_VIENNACL) && defined(PETSC_HAVE_VIENNACL_NO_CUDA)
5568:   PetscCall(MatSeqAIJRegister(MATMPIAIJVIENNACL, MatConvert_SeqAIJ_SeqAIJViennaCL));
5569: #endif
5570:   PetscFunctionReturn(PETSC_SUCCESS);
5571: }

5573: /*
5574:     Special version for direct calls from Fortran
5575: */
5576: #if defined(PETSC_HAVE_FORTRAN_CAPS)
5577:   #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
5578: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
5579:   #define matsetvaluesseqaij_ matsetvaluesseqaij
5580: #endif

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

5584: /* Change these macros so can be used in void function */
5585: /* Identical to PetscCallVoid, except it assigns to *_ierr */
5586: #undef PetscCall
5587: #define PetscCall(...) \
5588:   do { \
5589:     PetscErrorCode ierr_msv_mpiaij = __VA_ARGS__; \
5590:     if (PetscUnlikely(ierr_msv_mpiaij)) { \
5591:       *_ierr = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_msv_mpiaij, PETSC_ERROR_REPEAT, " "); \
5592:       return; \
5593:     } \
5594:   } while (0)

5596: #undef SETERRQ
5597: #define SETERRQ(comm, ierr, ...) \
5598:   do { \
5599:     *_ierr = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
5600:     return; \
5601:   } while (0)

5603: PETSC_EXTERN void matsetvaluesseqaij_(Mat *AA, PetscInt *mm, const PetscInt im[], PetscInt *nn, const PetscInt in[], const PetscScalar v[], InsertMode *isis, PetscErrorCode *_ierr)
5604: {
5605:   Mat         A = *AA;
5606:   PetscInt    m = *mm, n = *nn;
5607:   InsertMode  is = *isis;
5608:   Mat_SeqAIJ *a  = (Mat_SeqAIJ *)A->data;
5609:   PetscInt   *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N;
5610:   PetscInt   *imax, *ai, *ailen;
5611:   PetscInt   *aj, nonew = a->nonew, lastcol = -1;
5612:   MatScalar  *ap, value, *aa;
5613:   PetscBool   ignorezeroentries = a->ignorezeroentries;
5614:   PetscBool   roworiented       = a->roworiented;

5616:   PetscFunctionBegin;
5617:   MatCheckPreallocated(A, 1);
5618:   imax  = a->imax;
5619:   ai    = a->i;
5620:   ailen = a->ilen;
5621:   aj    = a->j;
5622:   aa    = a->a;

5624:   for (k = 0; k < m; k++) { /* loop over added rows */
5625:     row = im[k];
5626:     if (row < 0) continue;
5627:     PetscCheck(row < A->rmap->n, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Row too large");
5628:     rp   = aj + ai[row];
5629:     ap   = aa + ai[row];
5630:     rmax = imax[row];
5631:     nrow = ailen[row];
5632:     low  = 0;
5633:     high = nrow;
5634:     for (l = 0; l < n; l++) { /* loop over added columns */
5635:       if (in[l] < 0) continue;
5636:       PetscCheck(in[l] < A->cmap->n, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Column too large");
5637:       col = in[l];
5638:       if (roworiented) value = v[l + k * n];
5639:       else value = v[k + l * m];

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

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