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