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