Actual source code: diagonal.c
1: #include <petsc/private/matimpl.h>
3: typedef struct {
4: Vec diag;
5: PetscBool diag_valid;
6: Vec inv_diag;
7: PetscBool inv_diag_valid;
8: PetscObjectState diag_state, inv_diag_state;
9: PetscInt *col;
10: PetscScalar *val;
11: } Mat_Diagonal;
13: static PetscErrorCode MatDiagonalSetUpDiagonal(Mat A)
14: {
15: Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
17: PetscFunctionBegin;
18: if (!ctx->diag_valid) {
19: PetscAssert(ctx->inv_diag_valid, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Neither diagonal nor inverse diagonal is in a valid state");
20: PetscCall(VecCopy(ctx->inv_diag, ctx->diag));
21: PetscCall(VecReciprocal(ctx->diag));
22: ctx->diag_valid = PETSC_TRUE;
23: }
24: PetscFunctionReturn(PETSC_SUCCESS);
25: }
27: static PetscErrorCode MatDiagonalSetUpInverseDiagonal(Mat A)
28: {
29: Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
31: PetscFunctionBegin;
32: if (!ctx->inv_diag_valid) {
33: PetscAssert(ctx->diag_valid, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Neither diagonal nor inverse diagonal is in a valid state");
34: PetscCall(VecCopy(ctx->diag, ctx->inv_diag));
35: PetscCall(VecReciprocal(ctx->inv_diag));
36: ctx->inv_diag_valid = PETSC_TRUE;
37: }
38: PetscFunctionReturn(PETSC_SUCCESS);
39: }
41: static PetscErrorCode MatAXPY_Diagonal(Mat Y, PetscScalar a, Mat X, MatStructure str)
42: {
43: Mat_Diagonal *yctx = (Mat_Diagonal *)Y->data;
44: Mat_Diagonal *xctx = (Mat_Diagonal *)X->data;
46: PetscFunctionBegin;
47: PetscCall(MatDiagonalSetUpDiagonal(Y));
48: PetscCall(MatDiagonalSetUpDiagonal(X));
49: PetscCall(VecAXPY(yctx->diag, a, xctx->diag));
50: yctx->inv_diag_valid = PETSC_FALSE;
51: PetscFunctionReturn(PETSC_SUCCESS);
52: }
54: static PetscErrorCode MatGetRow_Diagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **vals)
55: {
56: Mat_Diagonal *mat = (Mat_Diagonal *)A->data;
58: PetscFunctionBegin;
59: if (ncols) *ncols = 1;
60: if (cols) {
61: if (!mat->col) PetscCall(PetscMalloc1(1, &mat->col));
62: *mat->col = row;
63: *cols = mat->col;
64: }
65: if (vals) {
66: const PetscScalar *v;
68: if (!mat->val) PetscCall(PetscMalloc1(1, &mat->val));
69: PetscCall(VecGetArrayRead(mat->diag, &v));
70: *mat->val = v[row];
71: *vals = mat->val;
72: PetscCall(VecRestoreArrayRead(mat->diag, &v));
73: }
74: PetscFunctionReturn(PETSC_SUCCESS);
75: }
77: static PetscErrorCode MatRestoreRow_Diagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **vals)
78: {
79: PetscFunctionBegin;
80: PetscFunctionReturn(PETSC_SUCCESS);
81: }
83: static PetscErrorCode MatMult_Diagonal(Mat A, Vec x, Vec y)
84: {
85: Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
87: PetscFunctionBegin;
88: PetscCall(MatDiagonalSetUpDiagonal(A));
89: PetscCall(VecPointwiseMult(y, ctx->diag, x));
90: PetscFunctionReturn(PETSC_SUCCESS);
91: }
93: static PetscErrorCode MatMultAdd_Diagonal(Mat mat, Vec v1, Vec v2, Vec v3)
94: {
95: Mat_Diagonal *ctx = (Mat_Diagonal *)mat->data;
97: PetscFunctionBegin;
98: PetscCall(MatDiagonalSetUpDiagonal(mat));
99: if (v2 != v3) {
100: PetscCall(VecPointwiseMult(v3, ctx->diag, v1));
101: PetscCall(VecAXPY(v3, 1.0, v2));
102: } else {
103: Vec w;
104: PetscCall(VecDuplicate(v3, &w));
105: PetscCall(VecPointwiseMult(w, ctx->diag, v1));
106: PetscCall(VecAXPY(v3, 1.0, w));
107: PetscCall(VecDestroy(&w));
108: }
109: PetscFunctionReturn(PETSC_SUCCESS);
110: }
112: static PetscErrorCode MatNorm_Diagonal(Mat A, NormType type, PetscReal *nrm)
113: {
114: Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
116: PetscFunctionBegin;
117: PetscCall(MatDiagonalSetUpDiagonal(A));
118: type = (type == NORM_FROBENIUS) ? NORM_2 : NORM_INFINITY;
119: PetscCall(VecNorm(ctx->diag, type, nrm));
120: PetscFunctionReturn(PETSC_SUCCESS);
121: }
123: static PetscErrorCode MatDuplicate_Diagonal(Mat A, MatDuplicateOption op, Mat *B)
124: {
125: Mat_Diagonal *actx = (Mat_Diagonal *)A->data;
127: PetscFunctionBegin;
128: PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
129: PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
130: PetscCall(MatSetBlockSizesFromMats(*B, A, A));
131: PetscCall(MatSetType(*B, MATDIAGONAL));
132: PetscCall(PetscLayoutReference(A->rmap, &(*B)->rmap));
133: PetscCall(PetscLayoutReference(A->cmap, &(*B)->cmap));
134: if (op == MAT_COPY_VALUES) {
135: Mat_Diagonal *bctx = (Mat_Diagonal *)(*B)->data;
137: PetscCall(MatSetUp(A));
138: PetscCall(MatSetUp(*B));
139: PetscCall(VecCopy(actx->diag, bctx->diag));
140: PetscCall(VecCopy(actx->inv_diag, bctx->inv_diag));
141: bctx->diag_valid = actx->diag_valid;
142: bctx->inv_diag_valid = actx->inv_diag_valid;
143: }
144: PetscFunctionReturn(PETSC_SUCCESS);
145: }
147: /*@
148: MatDiagonalGetDiagonal - Get the diagonal of a `MATDIAGONAL`
150: Input Parameter:
151: . A - the `MATDIAGONAL`
153: Output Parameter:
154: . diag - the `Vec` that defines the diagonal
156: Level: developer
158: Note:
159: The user must call
160: `MatDiagonalRestoreDiagonal()` before using the matrix again.
162: For a copy of the diagonal values, rather than a reference, use `MatGetDiagonal()`
164: Any changes to the obtained vector immediately change the action of the `Mat`. The matrix can be changed more efficiently by accessing this vector and changing its values, instead of filling a work vector and using `MatDiagonalSet()`
166: .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()`, `MatGetDiagonal()`
167: @*/
168: PetscErrorCode MatDiagonalGetDiagonal(Mat A, Vec *diag)
169: {
170: PetscFunctionBegin;
172: PetscAssertPointer(diag, 2);
173: *diag = NULL;
174: PetscUseMethod((PetscObject)A, "MatDiagonalGetDiagonal_C", (Mat, Vec *), (A, diag));
175: PetscFunctionReturn(PETSC_SUCCESS);
176: }
178: static PetscErrorCode MatDiagonalGetDiagonal_Diagonal(Mat A, Vec *diag)
179: {
180: Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
182: PetscFunctionBegin;
183: PetscCall(MatSetUp(A));
184: PetscCall(MatDiagonalSetUpDiagonal(A));
185: *diag = ctx->diag;
186: PetscCall(PetscObjectStateGet((PetscObject)*diag, &ctx->diag_state));
187: PetscFunctionReturn(PETSC_SUCCESS);
188: }
190: /*@
191: MatDiagonalRestoreDiagonal - Restore the diagonal of a `MATDIAGONAL`
193: Input Parameters:
194: + A - the `MATDIAGONAL`
195: - diag - the `Vec` obtained from `MatDiagonalGetDiagonal()`
197: Level: developer
199: Note:
200: Use `MatDiagonalSet()` to change the values by copy, rather than reference.
202: .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalGetDiagonal()`
203: @*/
204: PetscErrorCode MatDiagonalRestoreDiagonal(Mat A, Vec *diag)
205: {
206: PetscFunctionBegin;
208: PetscAssertPointer(diag, 2);
209: PetscUseMethod((PetscObject)A, "MatDiagonalRestoreDiagonal_C", (Mat, Vec *), (A, diag));
210: PetscFunctionReturn(PETSC_SUCCESS);
211: }
213: static PetscErrorCode MatDiagonalRestoreDiagonal_Diagonal(Mat A, Vec *diag)
214: {
215: Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
216: PetscObjectState diag_state;
218: PetscFunctionBegin;
219: PetscCheck(ctx->diag == *diag, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Restored a different diagonal vector");
220: ctx->diag_valid = PETSC_TRUE;
221: PetscCall(PetscObjectStateGet((PetscObject)*diag, &diag_state));
222: if (ctx->diag_state != diag_state) {
223: PetscCall(PetscObjectStateIncrease((PetscObject)A));
224: ctx->inv_diag_valid = PETSC_FALSE;
225: }
226: *diag = NULL;
227: PetscFunctionReturn(PETSC_SUCCESS);
228: }
230: /*@
231: MatDiagonalGetInverseDiagonal - Get the inverse diagonal of a `MATDIAGONAL`
233: Input Parameter:
234: . A - the `MATDIAGONAL`
236: Output Parameter:
237: . inv_diag - the `Vec` that defines the inverse diagonal
239: Level: developer
241: Note:
242: The user must call
243: `MatDiagonalRestoreInverseDiagonal()` before using the matrix again.
245: If a matrix is created only to call `MatSolve()` (which happens for `MATLMVMDIAGBROYDEN`),
246: using `MatDiagonalGetInverseDiagonal()` and `MatDiagonalRestoreInverseDiagonal()` avoids copies
247: and avoids any call to `VecReciprocal()`.
249: .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MATLMVMBROYDEN`, `MatSolve()`
250: @*/
251: PetscErrorCode MatDiagonalGetInverseDiagonal(Mat A, Vec *inv_diag)
252: {
253: PetscFunctionBegin;
255: PetscAssertPointer(inv_diag, 2);
256: *inv_diag = NULL;
257: PetscUseMethod((PetscObject)A, "MatDiagonalGetInverseDiagonal_C", (Mat, Vec *), (A, inv_diag));
258: PetscFunctionReturn(PETSC_SUCCESS);
259: }
261: static PetscErrorCode MatDiagonalGetInverseDiagonal_Diagonal(Mat A, Vec *inv_diag)
262: {
263: Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
265: PetscFunctionBegin;
266: PetscCall(MatSetUp(A));
267: PetscCall(MatDiagonalSetUpInverseDiagonal(A));
268: *inv_diag = ctx->inv_diag;
269: PetscCall(PetscObjectStateGet((PetscObject)*inv_diag, &ctx->inv_diag_state));
270: PetscFunctionReturn(PETSC_SUCCESS);
271: }
273: /*@
274: MatDiagonalRestoreInverseDiagonal - Restore the inverse diagonal of a `MATDIAGONAL`
276: Input Parameters:
277: + A - the `MATDIAGONAL`
278: - inv_diag - the `Vec` obtained from `MatDiagonalGetInverseDiagonal()`
280: Level: developer
282: .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalGetInverseDiagonal()`
283: @*/
284: PetscErrorCode MatDiagonalRestoreInverseDiagonal(Mat A, Vec *inv_diag)
285: {
286: PetscFunctionBegin;
288: PetscAssertPointer(inv_diag, 2);
289: PetscUseMethod((PetscObject)A, "MatDiagonalRestoreInverseDiagonal_C", (Mat, Vec *), (A, inv_diag));
290: PetscFunctionReturn(PETSC_SUCCESS);
291: }
293: static PetscErrorCode MatDiagonalRestoreInverseDiagonal_Diagonal(Mat A, Vec *inv_diag)
294: {
295: Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
296: PetscObjectState inv_diag_state;
298: PetscFunctionBegin;
299: PetscCheck(ctx->inv_diag == *inv_diag, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Restored a different diagonal vector");
300: ctx->inv_diag_valid = PETSC_TRUE;
301: PetscCall(PetscObjectStateGet((PetscObject)*inv_diag, &inv_diag_state));
302: if (ctx->inv_diag_state != inv_diag_state) {
303: PetscCall(PetscObjectStateIncrease((PetscObject)A));
304: ctx->diag_valid = PETSC_FALSE;
305: }
306: *inv_diag = NULL;
307: PetscFunctionReturn(PETSC_SUCCESS);
308: }
310: static PetscErrorCode MatDestroy_Diagonal(Mat mat)
311: {
312: Mat_Diagonal *ctx = (Mat_Diagonal *)mat->data;
314: PetscFunctionBegin;
315: PetscCall(VecDestroy(&ctx->diag));
316: PetscCall(VecDestroy(&ctx->inv_diag));
317: PetscCall(PetscFree(ctx->col));
318: PetscCall(PetscFree(ctx->val));
319: PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalGetDiagonal_C", NULL));
320: PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalRestoreDiagonal_C", NULL));
321: PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalGetInverseDiagonal_C", NULL));
322: PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalRestoreInverseDiagonal_C", NULL));
323: PetscCall(PetscFree(mat->data));
324: mat->structural_symmetry_eternal = PETSC_FALSE;
325: mat->symmetry_eternal = PETSC_FALSE;
326: PetscFunctionReturn(PETSC_SUCCESS);
327: }
329: static PetscErrorCode MatView_Diagonal(Mat J, PetscViewer viewer)
330: {
331: Mat_Diagonal *ctx = (Mat_Diagonal *)J->data;
332: PetscBool iascii;
334: PetscFunctionBegin;
335: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
336: if (iascii) {
337: PetscViewerFormat format;
339: PetscCall(PetscViewerGetFormat(viewer, &format));
340: if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS);
341: PetscCall(MatDiagonalSetUpDiagonal(J));
342: PetscCall(VecView(ctx->diag, viewer));
343: }
344: PetscFunctionReturn(PETSC_SUCCESS);
345: }
347: static PetscErrorCode MatGetDiagonal_Diagonal(Mat J, Vec x)
348: {
349: Mat_Diagonal *ctx = (Mat_Diagonal *)J->data;
351: PetscFunctionBegin;
352: PetscCall(MatDiagonalSetUpDiagonal(J));
353: PetscCall(VecCopy(ctx->diag, x));
354: PetscFunctionReturn(PETSC_SUCCESS);
355: }
357: static PetscErrorCode MatDiagonalSet_Diagonal(Mat J, Vec D, InsertMode is)
358: {
359: Mat_Diagonal *ctx = (Mat_Diagonal *)J->data;
361: PetscFunctionBegin;
362: switch (is) {
363: case ADD_VALUES:
364: case ADD_ALL_VALUES:
365: case ADD_BC_VALUES:
366: PetscCall(MatDiagonalSetUpDiagonal(J));
367: PetscCall(VecAXPY(ctx->diag, 1.0, D));
368: ctx->inv_diag_valid = PETSC_FALSE;
369: break;
370: case INSERT_VALUES:
371: case INSERT_BC_VALUES:
372: case INSERT_ALL_VALUES:
373: PetscCall(MatSetUp(J));
374: PetscCall(VecCopy(D, ctx->diag));
375: ctx->diag_valid = PETSC_TRUE;
376: ctx->inv_diag_valid = PETSC_FALSE;
377: break;
378: case MAX_VALUES:
379: PetscCall(MatDiagonalSetUpDiagonal(J));
380: PetscCall(VecPointwiseMax(ctx->diag, D, ctx->diag));
381: ctx->inv_diag_valid = PETSC_FALSE;
382: break;
383: case MIN_VALUES:
384: PetscCall(MatDiagonalSetUpDiagonal(J));
385: PetscCall(VecPointwiseMin(ctx->diag, D, ctx->diag));
386: ctx->inv_diag_valid = PETSC_FALSE;
387: break;
388: case NOT_SET_VALUES:
389: break;
390: }
391: PetscFunctionReturn(PETSC_SUCCESS);
392: }
394: static PetscErrorCode MatShift_Diagonal(Mat Y, PetscScalar a)
395: {
396: Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data;
398: PetscFunctionBegin;
399: PetscCall(MatDiagonalSetUpDiagonal(Y));
400: PetscCall(VecShift(ctx->diag, a));
401: ctx->inv_diag_valid = PETSC_FALSE;
402: PetscFunctionReturn(PETSC_SUCCESS);
403: }
405: static PetscErrorCode MatScale_Diagonal(Mat Y, PetscScalar a)
406: {
407: Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data;
409: PetscFunctionBegin;
410: PetscCall(MatSetUp(Y));
411: PetscCall(MatDiagonalSetUpDiagonal(Y));
412: PetscCall(VecScale(ctx->diag, a));
413: ctx->inv_diag_valid = PETSC_FALSE;
414: PetscFunctionReturn(PETSC_SUCCESS);
415: }
417: static PetscErrorCode MatDiagonalScale_Diagonal(Mat Y, Vec l, Vec r)
418: {
419: Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data;
421: PetscFunctionBegin;
422: PetscCall(MatDiagonalSetUpDiagonal(Y));
423: if (l) {
424: PetscCall(VecPointwiseMult(ctx->diag, ctx->diag, l));
425: ctx->inv_diag_valid = PETSC_FALSE;
426: }
427: if (r) {
428: PetscCall(VecPointwiseMult(ctx->diag, ctx->diag, r));
429: ctx->inv_diag_valid = PETSC_FALSE;
430: }
431: PetscFunctionReturn(PETSC_SUCCESS);
432: }
434: static PetscErrorCode MatSetUp_Diagonal(Mat A)
435: {
436: Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
438: PetscFunctionBegin;
439: if (!ctx->diag) {
440: PetscCall(PetscLayoutSetUp(A->cmap));
441: PetscCall(PetscLayoutSetUp(A->rmap));
442: PetscCall(MatCreateVecs(A, &ctx->diag, NULL));
443: PetscCall(VecDuplicate(ctx->diag, &ctx->inv_diag));
444: PetscCall(VecZeroEntries(ctx->diag));
445: ctx->diag_valid = PETSC_TRUE;
446: }
447: A->assembled = PETSC_TRUE;
448: PetscFunctionReturn(PETSC_SUCCESS);
449: }
451: static PetscErrorCode MatZeroEntries_Diagonal(Mat Y)
452: {
453: Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data;
455: PetscFunctionBegin;
456: PetscCall(MatSetUp(Y));
457: PetscCall(VecZeroEntries(ctx->diag));
458: ctx->diag_valid = PETSC_TRUE;
459: ctx->inv_diag_valid = PETSC_FALSE;
460: PetscFunctionReturn(PETSC_SUCCESS);
461: }
463: static PetscErrorCode MatSolve_Diagonal(Mat matin, Vec b, Vec x)
464: {
465: Mat_Diagonal *ctx = (Mat_Diagonal *)matin->data;
467: PetscFunctionBegin;
468: PetscCall(MatDiagonalSetUpInverseDiagonal(matin));
469: PetscCall(VecPointwiseMult(x, b, ctx->inv_diag));
470: PetscFunctionReturn(PETSC_SUCCESS);
471: }
473: static PetscErrorCode MatGetInfo_Diagonal(Mat A, MatInfoType flag, MatInfo *info)
474: {
475: PetscFunctionBegin;
476: info->block_size = 1.0;
477: info->nz_allocated = A->cmap->N;
478: info->nz_used = A->cmap->N;
479: info->nz_unneeded = 0.0;
480: info->assemblies = A->num_ass;
481: info->mallocs = 0.0;
482: info->memory = 0;
483: info->fill_ratio_given = 0;
484: info->fill_ratio_needed = 0;
485: info->factor_mallocs = 0;
486: PetscFunctionReturn(PETSC_SUCCESS);
487: }
489: /*@
490: MatCreateDiagonal - Creates a matrix defined by a given vector along its diagonal.
492: Collective
494: Input Parameter:
495: . diag - vector for the diagonal
497: Output Parameter:
498: . J - the diagonal matrix
500: Level: advanced
502: Notes:
503: Only supports square matrices with the same number of local rows and columns.
505: The input vector `diag` will be referenced internally: any changes to `diag`
506: will affect the matrix `J`.
508: .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatSolve()`
509: `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()`
510: @*/
511: PetscErrorCode MatCreateDiagonal(Vec diag, Mat *J)
512: {
513: PetscFunctionBegin;
515: PetscCall(MatCreate(PetscObjectComm((PetscObject)diag), J));
516: PetscInt m, M;
517: PetscCall(VecGetLocalSize(diag, &m));
518: PetscCall(VecGetSize(diag, &M));
519: PetscCall(MatSetSizes(*J, m, m, M, M));
520: PetscCall(MatSetType(*J, MATDIAGONAL));
522: PetscLayout map;
523: PetscCall(VecGetLayout(diag, &map));
524: PetscCall(MatSetLayouts(*J, map, map));
525: Mat_Diagonal *ctx = (Mat_Diagonal *)(*J)->data;
526: PetscCall(PetscObjectReference((PetscObject)diag));
527: PetscCall(VecDestroy(&ctx->diag));
528: PetscCall(VecDestroy(&ctx->inv_diag));
529: ctx->diag = diag;
530: ctx->diag_valid = PETSC_TRUE;
531: ctx->inv_diag_valid = PETSC_FALSE;
532: VecType type;
533: PetscCall(VecDuplicate(ctx->diag, &ctx->inv_diag));
534: PetscCall(VecGetType(diag, &type));
535: PetscCall(PetscFree((*J)->defaultvectype));
536: PetscCall(PetscStrallocpy(type, &(*J)->defaultvectype));
537: PetscCall(MatSetUp(*J));
538: ctx->col = NULL;
539: ctx->val = NULL;
540: PetscFunctionReturn(PETSC_SUCCESS);
541: }
543: /*MC
544: MATDIAGONAL - MATDIAGONAL = "diagonal" - A diagonal matrix type with the diagonal implemented as a `Vec`. Useful for
545: cases where `VecPointwiseMult()` or `VecPointwiseDivide()` should be thought of as the actions of a linear operator.
547: Level: advanced
549: .seealso: [](ch_matrices), `Mat`, `MatCreateDiagonal()`, `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()`
550: M*/
551: PETSC_INTERN PetscErrorCode MatCreate_Diagonal(Mat A)
552: {
553: Mat_Diagonal *ctx;
555: PetscFunctionBegin;
556: PetscCall(PetscNew(&ctx));
557: A->data = (void *)ctx;
559: A->structurally_symmetric = PETSC_BOOL3_TRUE;
560: A->structural_symmetry_eternal = PETSC_TRUE;
561: A->symmetry_eternal = PETSC_TRUE;
562: A->symmetric = PETSC_BOOL3_TRUE;
563: if (!PetscDefined(USE_COMPLEX)) A->hermitian = PETSC_BOOL3_TRUE;
565: A->ops->getrow = MatGetRow_Diagonal;
566: A->ops->restorerow = MatRestoreRow_Diagonal;
567: A->ops->mult = MatMult_Diagonal;
568: A->ops->multadd = MatMultAdd_Diagonal;
569: A->ops->multtranspose = MatMult_Diagonal;
570: A->ops->multtransposeadd = MatMultAdd_Diagonal;
571: A->ops->norm = MatNorm_Diagonal;
572: A->ops->duplicate = MatDuplicate_Diagonal;
573: A->ops->solve = MatSolve_Diagonal;
574: A->ops->solvetranspose = MatSolve_Diagonal;
575: A->ops->shift = MatShift_Diagonal;
576: A->ops->scale = MatScale_Diagonal;
577: A->ops->diagonalscale = MatDiagonalScale_Diagonal;
578: A->ops->getdiagonal = MatGetDiagonal_Diagonal;
579: A->ops->diagonalset = MatDiagonalSet_Diagonal;
580: A->ops->view = MatView_Diagonal;
581: A->ops->zeroentries = MatZeroEntries_Diagonal;
582: A->ops->destroy = MatDestroy_Diagonal;
583: A->ops->getinfo = MatGetInfo_Diagonal;
584: A->ops->axpy = MatAXPY_Diagonal;
585: A->ops->setup = MatSetUp_Diagonal;
587: PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalGetDiagonal_C", MatDiagonalGetDiagonal_Diagonal));
588: PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalRestoreDiagonal_C", MatDiagonalRestoreDiagonal_Diagonal));
589: PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalGetInverseDiagonal_C", MatDiagonalGetInverseDiagonal_Diagonal));
590: PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalRestoreInverseDiagonal_C", MatDiagonalRestoreInverseDiagonal_Diagonal));
591: PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATDIAGONAL));
592: PetscFunctionReturn(PETSC_SUCCESS);
593: }