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