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