xref: /petsc/src/mat/impls/diagonal/diagonal.c (revision ee8a0a41083610974418f8897c59118a81972f8d)
1345a4b08SToby Isaac #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
2345a4b08SToby Isaac 
3345a4b08SToby Isaac typedef struct {
4345a4b08SToby Isaac   Vec              diag;
5345a4b08SToby Isaac   PetscBool        diag_valid;
6345a4b08SToby Isaac   Vec              inv_diag;
7345a4b08SToby Isaac   PetscBool        inv_diag_valid;
8345a4b08SToby Isaac   PetscObjectState diag_state, inv_diag_state;
9c3e1b152SPierre Jolivet   PetscInt        *col;
10c3e1b152SPierre Jolivet   PetscScalar     *val;
11345a4b08SToby Isaac } Mat_Diagonal;
12345a4b08SToby Isaac 
13345a4b08SToby Isaac static PetscErrorCode MatDiagonalSetUpDiagonal(Mat A)
14345a4b08SToby Isaac {
15345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
16345a4b08SToby Isaac 
17345a4b08SToby Isaac   PetscFunctionBegin;
18345a4b08SToby Isaac   if (!ctx->diag_valid) {
19345a4b08SToby Isaac     PetscAssert(ctx->inv_diag_valid, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Neither diagonal nor inverse diagonal is in a valid state");
20345a4b08SToby Isaac     PetscCall(VecCopy(ctx->inv_diag, ctx->diag));
21345a4b08SToby Isaac     PetscCall(VecReciprocal(ctx->diag));
22345a4b08SToby Isaac     ctx->diag_valid = PETSC_TRUE;
23345a4b08SToby Isaac   }
24345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
25345a4b08SToby Isaac }
26345a4b08SToby Isaac 
27345a4b08SToby Isaac static PetscErrorCode MatDiagonalSetUpInverseDiagonal(Mat A)
28345a4b08SToby Isaac {
29345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
30345a4b08SToby Isaac 
31345a4b08SToby Isaac   PetscFunctionBegin;
32345a4b08SToby Isaac   if (!ctx->inv_diag_valid) {
33345a4b08SToby Isaac     PetscAssert(ctx->diag_valid, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Neither diagonal nor inverse diagonal is in a valid state");
34345a4b08SToby Isaac     PetscCall(VecCopy(ctx->diag, ctx->inv_diag));
35345a4b08SToby Isaac     PetscCall(VecReciprocal(ctx->inv_diag));
36345a4b08SToby Isaac     ctx->inv_diag_valid = PETSC_TRUE;
37345a4b08SToby Isaac   }
38345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
39345a4b08SToby Isaac }
40345a4b08SToby Isaac 
41345a4b08SToby Isaac static PetscErrorCode MatAXPY_Diagonal(Mat Y, PetscScalar a, Mat X, MatStructure str)
42345a4b08SToby Isaac {
43345a4b08SToby Isaac   Mat_Diagonal *yctx = (Mat_Diagonal *)Y->data;
44345a4b08SToby Isaac   Mat_Diagonal *xctx = (Mat_Diagonal *)X->data;
45345a4b08SToby Isaac 
46345a4b08SToby Isaac   PetscFunctionBegin;
47345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(Y));
48345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(X));
49345a4b08SToby Isaac   PetscCall(VecAXPY(yctx->diag, a, xctx->diag));
50345a4b08SToby Isaac   yctx->inv_diag_valid = PETSC_FALSE;
51345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
52345a4b08SToby Isaac }
53345a4b08SToby Isaac 
54c3e1b152SPierre Jolivet static PetscErrorCode MatGetRow_Diagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **vals)
55c3e1b152SPierre Jolivet {
56c3e1b152SPierre Jolivet   Mat_Diagonal *mat = (Mat_Diagonal *)A->data;
57c3e1b152SPierre Jolivet 
58c3e1b152SPierre Jolivet   PetscFunctionBegin;
59c3e1b152SPierre Jolivet   if (ncols) *ncols = 1;
60c3e1b152SPierre Jolivet   if (cols) {
61c3e1b152SPierre Jolivet     if (!mat->col) PetscCall(PetscMalloc1(1, &mat->col));
62c3e1b152SPierre Jolivet     *mat->col = row;
63c3e1b152SPierre Jolivet     *cols     = mat->col;
64c3e1b152SPierre Jolivet   }
65c3e1b152SPierre Jolivet   if (vals) {
66c3e1b152SPierre Jolivet     const PetscScalar *v;
67c3e1b152SPierre Jolivet 
68c3e1b152SPierre Jolivet     if (!mat->val) PetscCall(PetscMalloc1(1, &mat->val));
69c3e1b152SPierre Jolivet     PetscCall(VecGetArrayRead(mat->diag, &v));
70c3e1b152SPierre Jolivet     *mat->val = v[row];
71c3e1b152SPierre Jolivet     *vals     = mat->val;
72c3e1b152SPierre Jolivet     PetscCall(VecRestoreArrayRead(mat->diag, &v));
73c3e1b152SPierre Jolivet   }
74c3e1b152SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
75c3e1b152SPierre Jolivet }
76c3e1b152SPierre Jolivet 
77c3e1b152SPierre Jolivet static PetscErrorCode MatRestoreRow_Diagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **vals)
78c3e1b152SPierre Jolivet {
79c3e1b152SPierre Jolivet   PetscFunctionBegin;
80c3e1b152SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
81c3e1b152SPierre Jolivet }
82c3e1b152SPierre Jolivet 
83345a4b08SToby Isaac static PetscErrorCode MatMult_Diagonal(Mat A, Vec x, Vec y)
84345a4b08SToby Isaac {
85345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
86345a4b08SToby Isaac 
87345a4b08SToby Isaac   PetscFunctionBegin;
88345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(A));
89345a4b08SToby Isaac   PetscCall(VecPointwiseMult(y, ctx->diag, x));
90345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
91345a4b08SToby Isaac }
92345a4b08SToby Isaac 
93345a4b08SToby Isaac static PetscErrorCode MatMultAdd_Diagonal(Mat mat, Vec v1, Vec v2, Vec v3)
94345a4b08SToby Isaac {
95345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)mat->data;
96345a4b08SToby Isaac 
97345a4b08SToby Isaac   PetscFunctionBegin;
98345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(mat));
99345a4b08SToby Isaac   if (v2 != v3) {
100345a4b08SToby Isaac     PetscCall(VecPointwiseMult(v3, ctx->diag, v1));
101345a4b08SToby Isaac     PetscCall(VecAXPY(v3, 1.0, v2));
102345a4b08SToby Isaac   } else {
103345a4b08SToby Isaac     Vec w;
104345a4b08SToby Isaac     PetscCall(VecDuplicate(v3, &w));
105345a4b08SToby Isaac     PetscCall(VecPointwiseMult(w, ctx->diag, v1));
106345a4b08SToby Isaac     PetscCall(VecAXPY(v3, 1.0, w));
107345a4b08SToby Isaac     PetscCall(VecDestroy(&w));
108345a4b08SToby Isaac   }
109345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
110345a4b08SToby Isaac }
111345a4b08SToby Isaac 
112345a4b08SToby Isaac static PetscErrorCode MatNorm_Diagonal(Mat A, NormType type, PetscReal *nrm)
113345a4b08SToby Isaac {
114345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
115345a4b08SToby Isaac 
116345a4b08SToby Isaac   PetscFunctionBegin;
117345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(A));
118345a4b08SToby Isaac   type = (type == NORM_FROBENIUS) ? NORM_2 : NORM_INFINITY;
119345a4b08SToby Isaac   PetscCall(VecNorm(ctx->diag, type, nrm));
120345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
121345a4b08SToby Isaac }
122345a4b08SToby Isaac 
123345a4b08SToby Isaac static PetscErrorCode MatDuplicate_Diagonal(Mat A, MatDuplicateOption op, Mat *B)
124345a4b08SToby Isaac {
125345a4b08SToby Isaac   Mat_Diagonal *actx = (Mat_Diagonal *)A->data;
126345a4b08SToby Isaac 
127345a4b08SToby Isaac   PetscFunctionBegin;
128345a4b08SToby Isaac   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
129345a4b08SToby Isaac   PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
130345a4b08SToby Isaac   PetscCall(MatSetBlockSizesFromMats(*B, A, A));
131345a4b08SToby Isaac   PetscCall(MatSetType(*B, MATDIAGONAL));
132345a4b08SToby Isaac   PetscCall(PetscLayoutReference(A->rmap, &(*B)->rmap));
133345a4b08SToby Isaac   PetscCall(PetscLayoutReference(A->cmap, &(*B)->cmap));
134345a4b08SToby Isaac   if (op == MAT_COPY_VALUES) {
135345a4b08SToby Isaac     Mat_Diagonal *bctx = (Mat_Diagonal *)(*B)->data;
136345a4b08SToby Isaac 
137345a4b08SToby Isaac     PetscCall(MatSetUp(A));
138345a4b08SToby Isaac     PetscCall(MatSetUp(*B));
139345a4b08SToby Isaac     PetscCall(VecCopy(actx->diag, bctx->diag));
140345a4b08SToby Isaac     PetscCall(VecCopy(actx->inv_diag, bctx->inv_diag));
141345a4b08SToby Isaac     bctx->diag_valid     = actx->diag_valid;
142345a4b08SToby Isaac     bctx->inv_diag_valid = actx->inv_diag_valid;
143345a4b08SToby Isaac   }
144345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
145345a4b08SToby Isaac }
146345a4b08SToby Isaac 
147345a4b08SToby Isaac /*@
148345a4b08SToby Isaac   MatDiagonalGetDiagonal - Get the diagonal of a `MATDIAGONAL`
149345a4b08SToby Isaac 
150345a4b08SToby Isaac   Input Parameter:
151345a4b08SToby Isaac . A - the `MATDIAGONAL`
152345a4b08SToby Isaac 
153345a4b08SToby Isaac   Output Parameter:
154345a4b08SToby Isaac . diag - the `Vec` that defines the diagonal
155345a4b08SToby Isaac 
156345a4b08SToby Isaac   Level: developer
157345a4b08SToby Isaac 
158345a4b08SToby Isaac   Note:
159345a4b08SToby Isaac   The user must call
160345a4b08SToby Isaac   `MatDiagonalRestoreDiagonal()` before using the matrix again.
161345a4b08SToby Isaac 
162345a4b08SToby Isaac   For a copy of the diagonal values, rather than a reference, use `MatGetDiagonal()`
163345a4b08SToby Isaac 
164345a4b08SToby Isaac   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()`
165345a4b08SToby Isaac 
166be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()`, `MatGetDiagonal()`
167345a4b08SToby Isaac @*/
168345a4b08SToby Isaac PetscErrorCode MatDiagonalGetDiagonal(Mat A, Vec *diag)
169345a4b08SToby Isaac {
170345a4b08SToby Isaac   PetscFunctionBegin;
171345a4b08SToby Isaac   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1724f572ea9SToby Isaac   PetscAssertPointer(diag, 2);
173345a4b08SToby Isaac   *diag = NULL;
174345a4b08SToby Isaac   PetscUseMethod((PetscObject)A, "MatDiagonalGetDiagonal_C", (Mat, Vec *), (A, diag));
175345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
176345a4b08SToby Isaac }
177345a4b08SToby Isaac 
178345a4b08SToby Isaac static PetscErrorCode MatDiagonalGetDiagonal_Diagonal(Mat A, Vec *diag)
179345a4b08SToby Isaac {
180345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
181345a4b08SToby Isaac 
182345a4b08SToby Isaac   PetscFunctionBegin;
183345a4b08SToby Isaac   PetscCall(MatSetUp(A));
184345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(A));
185345a4b08SToby Isaac   *diag = ctx->diag;
186345a4b08SToby Isaac   PetscCall(PetscObjectStateGet((PetscObject)*diag, &ctx->diag_state));
187345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
188345a4b08SToby Isaac }
189345a4b08SToby Isaac 
190345a4b08SToby Isaac /*@
191345a4b08SToby Isaac   MatDiagonalRestoreDiagonal - Restore the diagonal of a `MATDIAGONAL`
192345a4b08SToby Isaac 
193345a4b08SToby Isaac   Input Parameters:
194345a4b08SToby Isaac + A    - the `MATDIAGONAL`
195345a4b08SToby Isaac - diag - the `Vec` obtained from `MatDiagonalGetDiagonal()`
196345a4b08SToby Isaac 
197345a4b08SToby Isaac   Level: developer
198345a4b08SToby Isaac 
199345a4b08SToby Isaac   Note:
200345a4b08SToby Isaac   Use `MatDiagonalSet()` to change the values by copy, rather than reference.
201345a4b08SToby Isaac 
202be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalGetDiagonal()`
203345a4b08SToby Isaac @*/
204345a4b08SToby Isaac PetscErrorCode MatDiagonalRestoreDiagonal(Mat A, Vec *diag)
205345a4b08SToby Isaac {
206345a4b08SToby Isaac   PetscFunctionBegin;
207345a4b08SToby Isaac   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
2084f572ea9SToby Isaac   PetscAssertPointer(diag, 2);
209345a4b08SToby Isaac   PetscUseMethod((PetscObject)A, "MatDiagonalRestoreDiagonal_C", (Mat, Vec *), (A, diag));
210345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
211345a4b08SToby Isaac }
212345a4b08SToby Isaac 
213345a4b08SToby Isaac static PetscErrorCode MatDiagonalRestoreDiagonal_Diagonal(Mat A, Vec *diag)
214345a4b08SToby Isaac {
215345a4b08SToby Isaac   Mat_Diagonal    *ctx = (Mat_Diagonal *)A->data;
216345a4b08SToby Isaac   PetscObjectState diag_state;
217345a4b08SToby Isaac 
218345a4b08SToby Isaac   PetscFunctionBegin;
219345a4b08SToby Isaac   PetscCheck(ctx->diag == *diag, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Restored a different diagonal vector");
220345a4b08SToby Isaac   ctx->diag_valid = PETSC_TRUE;
221345a4b08SToby Isaac   PetscCall(PetscObjectStateGet((PetscObject)*diag, &diag_state));
222345a4b08SToby Isaac   if (ctx->diag_state != diag_state) {
223345a4b08SToby Isaac     PetscCall(PetscObjectStateIncrease((PetscObject)A));
224345a4b08SToby Isaac     ctx->inv_diag_valid = PETSC_FALSE;
225345a4b08SToby Isaac   }
226345a4b08SToby Isaac   *diag = NULL;
227345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
228345a4b08SToby Isaac }
229345a4b08SToby Isaac 
230345a4b08SToby Isaac /*@
231345a4b08SToby Isaac   MatDiagonalGetInverseDiagonal - Get the inverse diagonal of a `MATDIAGONAL`
232345a4b08SToby Isaac 
233345a4b08SToby Isaac   Input Parameter:
234345a4b08SToby Isaac . A - the `MATDIAGONAL`
235345a4b08SToby Isaac 
236345a4b08SToby Isaac   Output Parameter:
237345a4b08SToby Isaac . inv_diag - the `Vec` that defines the inverse diagonal
238345a4b08SToby Isaac 
239345a4b08SToby Isaac   Level: developer
240345a4b08SToby Isaac 
241345a4b08SToby Isaac   Note:
242345a4b08SToby Isaac   The user must call
243345a4b08SToby Isaac   `MatDiagonalRestoreInverseDiagonal()` before using the matrix again.
244345a4b08SToby Isaac 
245345a4b08SToby Isaac   If a matrix is created only to call `MatSolve()` (which happens for `MATLMVMDIAGBROYDEN`),
246345a4b08SToby Isaac   using `MatDiagonalGetInverseDiagonal()` and `MatDiagonalRestoreInverseDiagonal()` avoids copies
247345a4b08SToby Isaac   and avoids any call to `VecReciprocal()`.
248345a4b08SToby Isaac 
249be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MATLMVMBROYDEN`, `MatSolve()`
250345a4b08SToby Isaac @*/
251345a4b08SToby Isaac PetscErrorCode MatDiagonalGetInverseDiagonal(Mat A, Vec *inv_diag)
252345a4b08SToby Isaac {
253345a4b08SToby Isaac   PetscFunctionBegin;
254345a4b08SToby Isaac   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
2554f572ea9SToby Isaac   PetscAssertPointer(inv_diag, 2);
256345a4b08SToby Isaac   *inv_diag = NULL;
257345a4b08SToby Isaac   PetscUseMethod((PetscObject)A, "MatDiagonalGetInverseDiagonal_C", (Mat, Vec *), (A, inv_diag));
258345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
259345a4b08SToby Isaac }
260345a4b08SToby Isaac 
261345a4b08SToby Isaac static PetscErrorCode MatDiagonalGetInverseDiagonal_Diagonal(Mat A, Vec *inv_diag)
262345a4b08SToby Isaac {
263345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
264345a4b08SToby Isaac 
265345a4b08SToby Isaac   PetscFunctionBegin;
266345a4b08SToby Isaac   PetscCall(MatSetUp(A));
267345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpInverseDiagonal(A));
268345a4b08SToby Isaac   *inv_diag = ctx->inv_diag;
269345a4b08SToby Isaac   PetscCall(PetscObjectStateGet((PetscObject)*inv_diag, &ctx->inv_diag_state));
270345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
271345a4b08SToby Isaac }
272345a4b08SToby Isaac 
273345a4b08SToby Isaac /*@
274345a4b08SToby Isaac   MatDiagonalRestoreInverseDiagonal - Restore the inverse diagonal of a `MATDIAGONAL`
275345a4b08SToby Isaac 
276345a4b08SToby Isaac   Input Parameters:
277345a4b08SToby Isaac + A        - the `MATDIAGONAL`
278345a4b08SToby Isaac - inv_diag - the `Vec` obtained from `MatDiagonalGetInverseDiagonal()`
279345a4b08SToby Isaac 
280345a4b08SToby Isaac   Level: developer
281345a4b08SToby Isaac 
282be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalGetInverseDiagonal()`
283345a4b08SToby Isaac @*/
284345a4b08SToby Isaac PetscErrorCode MatDiagonalRestoreInverseDiagonal(Mat A, Vec *inv_diag)
285345a4b08SToby Isaac {
286345a4b08SToby Isaac   PetscFunctionBegin;
287345a4b08SToby Isaac   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
2884f572ea9SToby Isaac   PetscAssertPointer(inv_diag, 2);
289345a4b08SToby Isaac   PetscUseMethod((PetscObject)A, "MatDiagonalRestoreInverseDiagonal_C", (Mat, Vec *), (A, inv_diag));
290345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
291345a4b08SToby Isaac }
292345a4b08SToby Isaac 
293345a4b08SToby Isaac static PetscErrorCode MatDiagonalRestoreInverseDiagonal_Diagonal(Mat A, Vec *inv_diag)
294345a4b08SToby Isaac {
295345a4b08SToby Isaac   Mat_Diagonal    *ctx = (Mat_Diagonal *)A->data;
296345a4b08SToby Isaac   PetscObjectState inv_diag_state;
297345a4b08SToby Isaac 
298345a4b08SToby Isaac   PetscFunctionBegin;
299345a4b08SToby Isaac   PetscCheck(ctx->inv_diag == *inv_diag, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Restored a different diagonal vector");
300345a4b08SToby Isaac   ctx->inv_diag_valid = PETSC_TRUE;
301345a4b08SToby Isaac   PetscCall(PetscObjectStateGet((PetscObject)*inv_diag, &inv_diag_state));
302345a4b08SToby Isaac   if (ctx->inv_diag_state != inv_diag_state) {
303345a4b08SToby Isaac     PetscCall(PetscObjectStateIncrease((PetscObject)A));
304345a4b08SToby Isaac     ctx->diag_valid = PETSC_FALSE;
305345a4b08SToby Isaac   }
306345a4b08SToby Isaac   *inv_diag = NULL;
307345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
308345a4b08SToby Isaac }
309345a4b08SToby Isaac 
310*ee8a0a41SPierre Jolivet static PetscErrorCode MatPermute_Diagonal(Mat A, IS rowp, IS colp, Mat *B)
311*ee8a0a41SPierre Jolivet {
312*ee8a0a41SPierre Jolivet   Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
313*ee8a0a41SPierre Jolivet   Vec           v;
314*ee8a0a41SPierre Jolivet 
315*ee8a0a41SPierre Jolivet   PetscFunctionBegin;
316*ee8a0a41SPierre Jolivet   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
317*ee8a0a41SPierre Jolivet   PetscCall(VecDuplicate(ctx->diag, &v));
318*ee8a0a41SPierre Jolivet   PetscCall(VecCopy(ctx->diag, v));
319*ee8a0a41SPierre Jolivet   PetscCall(VecPermute(v, rowp, PETSC_FALSE));
320*ee8a0a41SPierre Jolivet   PetscCall(MatCreateDiagonal(v, B));
321*ee8a0a41SPierre Jolivet   PetscCall(VecDestroy(&v));
322*ee8a0a41SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
323*ee8a0a41SPierre Jolivet }
324*ee8a0a41SPierre Jolivet 
325345a4b08SToby Isaac static PetscErrorCode MatDestroy_Diagonal(Mat mat)
326345a4b08SToby Isaac {
327345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)mat->data;
328345a4b08SToby Isaac 
329345a4b08SToby Isaac   PetscFunctionBegin;
330345a4b08SToby Isaac   PetscCall(VecDestroy(&ctx->diag));
331345a4b08SToby Isaac   PetscCall(VecDestroy(&ctx->inv_diag));
332c3e1b152SPierre Jolivet   PetscCall(PetscFree(ctx->col));
333c3e1b152SPierre Jolivet   PetscCall(PetscFree(ctx->val));
334345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalGetDiagonal_C", NULL));
335345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalRestoreDiagonal_C", NULL));
336345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalGetInverseDiagonal_C", NULL));
337345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalRestoreInverseDiagonal_C", NULL));
338e27ab85cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_diagonal_seqdense_C", NULL));
339e27ab85cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_diagonal_mpidense_C", NULL));
340345a4b08SToby Isaac   PetscCall(PetscFree(mat->data));
341345a4b08SToby Isaac   mat->structural_symmetry_eternal = PETSC_FALSE;
342345a4b08SToby Isaac   mat->symmetry_eternal            = PETSC_FALSE;
343345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
344345a4b08SToby Isaac }
345345a4b08SToby Isaac 
346345a4b08SToby Isaac static PetscErrorCode MatView_Diagonal(Mat J, PetscViewer viewer)
347345a4b08SToby Isaac {
348345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)J->data;
349345a4b08SToby Isaac   PetscBool     iascii;
350345a4b08SToby Isaac 
351345a4b08SToby Isaac   PetscFunctionBegin;
352345a4b08SToby Isaac   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
353345a4b08SToby Isaac   if (iascii) {
354345a4b08SToby Isaac     PetscViewerFormat format;
355345a4b08SToby Isaac 
356345a4b08SToby Isaac     PetscCall(PetscViewerGetFormat(viewer, &format));
357345a4b08SToby Isaac     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS);
358345a4b08SToby Isaac     PetscCall(MatDiagonalSetUpDiagonal(J));
359345a4b08SToby Isaac     PetscCall(VecView(ctx->diag, viewer));
360345a4b08SToby Isaac   }
361345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
362345a4b08SToby Isaac }
363345a4b08SToby Isaac 
364345a4b08SToby Isaac static PetscErrorCode MatGetDiagonal_Diagonal(Mat J, Vec x)
365345a4b08SToby Isaac {
366345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)J->data;
367345a4b08SToby Isaac 
368345a4b08SToby Isaac   PetscFunctionBegin;
369345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(J));
370345a4b08SToby Isaac   PetscCall(VecCopy(ctx->diag, x));
371345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
372345a4b08SToby Isaac }
373345a4b08SToby Isaac 
374345a4b08SToby Isaac static PetscErrorCode MatDiagonalSet_Diagonal(Mat J, Vec D, InsertMode is)
375345a4b08SToby Isaac {
376345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)J->data;
377345a4b08SToby Isaac 
378345a4b08SToby Isaac   PetscFunctionBegin;
379345a4b08SToby Isaac   switch (is) {
380345a4b08SToby Isaac   case ADD_VALUES:
381345a4b08SToby Isaac   case ADD_ALL_VALUES:
382345a4b08SToby Isaac   case ADD_BC_VALUES:
383345a4b08SToby Isaac     PetscCall(MatDiagonalSetUpDiagonal(J));
384345a4b08SToby Isaac     PetscCall(VecAXPY(ctx->diag, 1.0, D));
385345a4b08SToby Isaac     ctx->inv_diag_valid = PETSC_FALSE;
386345a4b08SToby Isaac     break;
387345a4b08SToby Isaac   case INSERT_VALUES:
388345a4b08SToby Isaac   case INSERT_BC_VALUES:
389345a4b08SToby Isaac   case INSERT_ALL_VALUES:
390345a4b08SToby Isaac     PetscCall(MatSetUp(J));
391345a4b08SToby Isaac     PetscCall(VecCopy(D, ctx->diag));
392345a4b08SToby Isaac     ctx->diag_valid     = PETSC_TRUE;
393345a4b08SToby Isaac     ctx->inv_diag_valid = PETSC_FALSE;
394345a4b08SToby Isaac     break;
395345a4b08SToby Isaac   case MAX_VALUES:
396345a4b08SToby Isaac     PetscCall(MatDiagonalSetUpDiagonal(J));
397345a4b08SToby Isaac     PetscCall(VecPointwiseMax(ctx->diag, D, ctx->diag));
398345a4b08SToby Isaac     ctx->inv_diag_valid = PETSC_FALSE;
399345a4b08SToby Isaac     break;
400345a4b08SToby Isaac   case MIN_VALUES:
401345a4b08SToby Isaac     PetscCall(MatDiagonalSetUpDiagonal(J));
402345a4b08SToby Isaac     PetscCall(VecPointwiseMin(ctx->diag, D, ctx->diag));
403345a4b08SToby Isaac     ctx->inv_diag_valid = PETSC_FALSE;
404345a4b08SToby Isaac     break;
405345a4b08SToby Isaac   case NOT_SET_VALUES:
406345a4b08SToby Isaac     break;
407345a4b08SToby Isaac   }
408345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
409345a4b08SToby Isaac }
410345a4b08SToby Isaac 
411345a4b08SToby Isaac static PetscErrorCode MatShift_Diagonal(Mat Y, PetscScalar a)
412345a4b08SToby Isaac {
413345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data;
414345a4b08SToby Isaac 
415345a4b08SToby Isaac   PetscFunctionBegin;
416345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(Y));
417345a4b08SToby Isaac   PetscCall(VecShift(ctx->diag, a));
418345a4b08SToby Isaac   ctx->inv_diag_valid = PETSC_FALSE;
419345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
420345a4b08SToby Isaac }
421345a4b08SToby Isaac 
422345a4b08SToby Isaac static PetscErrorCode MatScale_Diagonal(Mat Y, PetscScalar a)
423345a4b08SToby Isaac {
424345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data;
425345a4b08SToby Isaac 
426345a4b08SToby Isaac   PetscFunctionBegin;
427345a4b08SToby Isaac   PetscCall(MatSetUp(Y));
428345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(Y));
429345a4b08SToby Isaac   PetscCall(VecScale(ctx->diag, a));
430345a4b08SToby Isaac   ctx->inv_diag_valid = PETSC_FALSE;
431345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
432345a4b08SToby Isaac }
433345a4b08SToby Isaac 
434345a4b08SToby Isaac static PetscErrorCode MatDiagonalScale_Diagonal(Mat Y, Vec l, Vec r)
435345a4b08SToby Isaac {
436345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data;
437345a4b08SToby Isaac 
438345a4b08SToby Isaac   PetscFunctionBegin;
439345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpDiagonal(Y));
440345a4b08SToby Isaac   if (l) {
441345a4b08SToby Isaac     PetscCall(VecPointwiseMult(ctx->diag, ctx->diag, l));
442345a4b08SToby Isaac     ctx->inv_diag_valid = PETSC_FALSE;
443345a4b08SToby Isaac   }
444345a4b08SToby Isaac   if (r) {
445345a4b08SToby Isaac     PetscCall(VecPointwiseMult(ctx->diag, ctx->diag, r));
446345a4b08SToby Isaac     ctx->inv_diag_valid = PETSC_FALSE;
447345a4b08SToby Isaac   }
448345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
449345a4b08SToby Isaac }
450345a4b08SToby Isaac 
451345a4b08SToby Isaac static PetscErrorCode MatSetUp_Diagonal(Mat A)
452345a4b08SToby Isaac {
453345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)A->data;
454345a4b08SToby Isaac 
455345a4b08SToby Isaac   PetscFunctionBegin;
456345a4b08SToby Isaac   if (!ctx->diag) {
457345a4b08SToby Isaac     PetscCall(PetscLayoutSetUp(A->cmap));
458345a4b08SToby Isaac     PetscCall(PetscLayoutSetUp(A->rmap));
459345a4b08SToby Isaac     PetscCall(MatCreateVecs(A, &ctx->diag, NULL));
460345a4b08SToby Isaac     PetscCall(VecDuplicate(ctx->diag, &ctx->inv_diag));
461345a4b08SToby Isaac     PetscCall(VecZeroEntries(ctx->diag));
462345a4b08SToby Isaac     ctx->diag_valid = PETSC_TRUE;
463345a4b08SToby Isaac   }
464345a4b08SToby Isaac   A->assembled = PETSC_TRUE;
465345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
466345a4b08SToby Isaac }
467345a4b08SToby Isaac 
468345a4b08SToby Isaac static PetscErrorCode MatZeroEntries_Diagonal(Mat Y)
469345a4b08SToby Isaac {
470345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data;
471345a4b08SToby Isaac 
472345a4b08SToby Isaac   PetscFunctionBegin;
473345a4b08SToby Isaac   PetscCall(MatSetUp(Y));
474345a4b08SToby Isaac   PetscCall(VecZeroEntries(ctx->diag));
475345a4b08SToby Isaac   ctx->diag_valid     = PETSC_TRUE;
476345a4b08SToby Isaac   ctx->inv_diag_valid = PETSC_FALSE;
477345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
478345a4b08SToby Isaac }
479345a4b08SToby Isaac 
48066976f2fSJacob Faibussowitsch static PetscErrorCode MatSolve_Diagonal(Mat matin, Vec b, Vec x)
481345a4b08SToby Isaac {
482345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)matin->data;
483345a4b08SToby Isaac 
484345a4b08SToby Isaac   PetscFunctionBegin;
485345a4b08SToby Isaac   PetscCall(MatDiagonalSetUpInverseDiagonal(matin));
486345a4b08SToby Isaac   PetscCall(VecPointwiseMult(x, b, ctx->inv_diag));
487345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
488345a4b08SToby Isaac }
489345a4b08SToby Isaac 
49066976f2fSJacob Faibussowitsch static PetscErrorCode MatGetInfo_Diagonal(Mat A, MatInfoType flag, MatInfo *info)
491345a4b08SToby Isaac {
492345a4b08SToby Isaac   PetscFunctionBegin;
493345a4b08SToby Isaac   info->block_size        = 1.0;
494345a4b08SToby Isaac   info->nz_allocated      = A->cmap->N;
495345a4b08SToby Isaac   info->nz_used           = A->cmap->N;
496345a4b08SToby Isaac   info->nz_unneeded       = 0.0;
497345a4b08SToby Isaac   info->assemblies        = A->num_ass;
498345a4b08SToby Isaac   info->mallocs           = 0.0;
499345a4b08SToby Isaac   info->memory            = 0;
500345a4b08SToby Isaac   info->fill_ratio_given  = 0;
501345a4b08SToby Isaac   info->fill_ratio_needed = 0;
502345a4b08SToby Isaac   info->factor_mallocs    = 0;
503345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
504345a4b08SToby Isaac }
505345a4b08SToby Isaac 
506345a4b08SToby Isaac /*@
507345a4b08SToby Isaac   MatCreateDiagonal - Creates a matrix defined by a given vector along its diagonal.
508345a4b08SToby Isaac 
509345a4b08SToby Isaac   Collective
510345a4b08SToby Isaac 
511345a4b08SToby Isaac   Input Parameter:
512345a4b08SToby Isaac . diag - vector for the diagonal
513345a4b08SToby Isaac 
514345a4b08SToby Isaac   Output Parameter:
515345a4b08SToby Isaac . J - the diagonal matrix
516345a4b08SToby Isaac 
517345a4b08SToby Isaac   Level: advanced
518345a4b08SToby Isaac 
519345a4b08SToby Isaac   Notes:
520345a4b08SToby Isaac   Only supports square matrices with the same number of local rows and columns.
521345a4b08SToby Isaac 
522345a4b08SToby Isaac   The input vector `diag` will be referenced internally: any changes to `diag`
523345a4b08SToby Isaac   will affect the matrix `J`.
524345a4b08SToby Isaac 
525be50c303SSatish Balay .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatSolve()`
526345a4b08SToby Isaac           `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()`
527345a4b08SToby Isaac @*/
528345a4b08SToby Isaac PetscErrorCode MatCreateDiagonal(Vec diag, Mat *J)
529345a4b08SToby Isaac {
530345a4b08SToby Isaac   PetscFunctionBegin;
531345a4b08SToby Isaac   PetscValidHeaderSpecific(diag, VEC_CLASSID, 1);
532345a4b08SToby Isaac   PetscCall(MatCreate(PetscObjectComm((PetscObject)diag), J));
533345a4b08SToby Isaac   PetscInt m, M;
534345a4b08SToby Isaac   PetscCall(VecGetLocalSize(diag, &m));
535345a4b08SToby Isaac   PetscCall(VecGetSize(diag, &M));
536345a4b08SToby Isaac   PetscCall(MatSetSizes(*J, m, m, M, M));
537345a4b08SToby Isaac   PetscCall(MatSetType(*J, MATDIAGONAL));
538345a4b08SToby Isaac 
539345a4b08SToby Isaac   PetscLayout map;
540345a4b08SToby Isaac   PetscCall(VecGetLayout(diag, &map));
541345a4b08SToby Isaac   PetscCall(MatSetLayouts(*J, map, map));
542345a4b08SToby Isaac   Mat_Diagonal *ctx = (Mat_Diagonal *)(*J)->data;
543345a4b08SToby Isaac   PetscCall(PetscObjectReference((PetscObject)diag));
544345a4b08SToby Isaac   PetscCall(VecDestroy(&ctx->diag));
545345a4b08SToby Isaac   PetscCall(VecDestroy(&ctx->inv_diag));
546345a4b08SToby Isaac   ctx->diag           = diag;
547345a4b08SToby Isaac   ctx->diag_valid     = PETSC_TRUE;
548345a4b08SToby Isaac   ctx->inv_diag_valid = PETSC_FALSE;
549345a4b08SToby Isaac   VecType type;
550345a4b08SToby Isaac   PetscCall(VecDuplicate(ctx->diag, &ctx->inv_diag));
551345a4b08SToby Isaac   PetscCall(VecGetType(diag, &type));
552345a4b08SToby Isaac   PetscCall(PetscFree((*J)->defaultvectype));
553345a4b08SToby Isaac   PetscCall(PetscStrallocpy(type, &(*J)->defaultvectype));
554345a4b08SToby Isaac   PetscCall(MatSetUp(*J));
555c3e1b152SPierre Jolivet   ctx->col = NULL;
556c3e1b152SPierre Jolivet   ctx->val = NULL;
557345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
558345a4b08SToby Isaac }
559345a4b08SToby Isaac 
560e27ab85cSPierre Jolivet static PetscErrorCode MatProductNumeric_Diagonal_Dense(Mat C)
561e27ab85cSPierre Jolivet {
562e27ab85cSPierre Jolivet   Mat                A, B;
563e27ab85cSPierre Jolivet   Mat_Diagonal      *a;
564e27ab85cSPierre Jolivet   PetscScalar       *c;
565e27ab85cSPierre Jolivet   const PetscScalar *b, *alpha;
566e27ab85cSPierre Jolivet   PetscInt           ldb, ldc;
567e27ab85cSPierre Jolivet 
568e27ab85cSPierre Jolivet   PetscFunctionBegin;
569e27ab85cSPierre Jolivet   MatCheckProduct(C, 1);
570e27ab85cSPierre Jolivet   A = C->product->A;
571e27ab85cSPierre Jolivet   B = C->product->B;
572e27ab85cSPierre Jolivet   a = (Mat_Diagonal *)A->data;
573e27ab85cSPierre Jolivet   PetscCall(VecGetArrayRead(a->diag, &alpha));
574e27ab85cSPierre Jolivet   PetscCall(MatDenseGetLDA(B, &ldb));
575e27ab85cSPierre Jolivet   PetscCall(MatDenseGetLDA(C, &ldc));
576e27ab85cSPierre Jolivet   PetscCall(MatDenseGetArrayRead(B, &b));
577e27ab85cSPierre Jolivet   PetscCall(MatDenseGetArrayWrite(C, &c));
578e27ab85cSPierre Jolivet   for (PetscInt j = 0; j < B->cmap->N; j++)
579e27ab85cSPierre Jolivet     for (PetscInt i = 0; i < B->rmap->n; i++) c[i + j * ldc] = alpha[i] * b[i + j * ldb];
580e27ab85cSPierre Jolivet   PetscCall(MatDenseRestoreArrayWrite(C, &c));
581e27ab85cSPierre Jolivet   PetscCall(MatDenseRestoreArrayRead(B, &b));
582e27ab85cSPierre Jolivet   PetscCall(VecRestoreArrayRead(a->diag, &alpha));
583e27ab85cSPierre Jolivet   PetscCall(PetscLogFlops(B->cmap->N * B->rmap->n));
584e27ab85cSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
585e27ab85cSPierre Jolivet }
586e27ab85cSPierre Jolivet 
587e27ab85cSPierre Jolivet static PetscErrorCode MatProductSymbolic_Diagonal_Dense(Mat C)
588e27ab85cSPierre Jolivet {
589e27ab85cSPierre Jolivet   Mat      A, B;
590e27ab85cSPierre Jolivet   PetscInt n, N, m, M;
591e27ab85cSPierre Jolivet 
592e27ab85cSPierre Jolivet   PetscFunctionBegin;
593e27ab85cSPierre Jolivet   MatCheckProduct(C, 1);
594e27ab85cSPierre Jolivet   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
595e27ab85cSPierre Jolivet   A = C->product->A;
596e27ab85cSPierre Jolivet   B = C->product->B;
597e27ab85cSPierre Jolivet   PetscCall(MatDiagonalSetUpDiagonal(A));
598e27ab85cSPierre Jolivet   PetscCall(MatGetLocalSize(C, &m, &n));
599e27ab85cSPierre Jolivet   PetscCall(MatGetSize(C, &M, &N));
600e27ab85cSPierre Jolivet   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
601e27ab85cSPierre Jolivet     PetscCall(MatGetLocalSize(B, NULL, &n));
602e27ab85cSPierre Jolivet     PetscCall(MatGetSize(B, NULL, &N));
603e27ab85cSPierre Jolivet     PetscCall(MatGetLocalSize(A, &m, NULL));
604e27ab85cSPierre Jolivet     PetscCall(MatGetSize(A, &M, NULL));
605e27ab85cSPierre Jolivet     PetscCall(MatSetSizes(C, m, n, M, N));
606e27ab85cSPierre Jolivet   }
607e27ab85cSPierre Jolivet   PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
608e27ab85cSPierre Jolivet   PetscCall(MatSetUp(C));
609e27ab85cSPierre Jolivet   C->ops->productnumeric = MatProductNumeric_Diagonal_Dense;
610e27ab85cSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
611e27ab85cSPierre Jolivet }
612e27ab85cSPierre Jolivet 
613e27ab85cSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Diagonal_Dense_AB(Mat C)
614e27ab85cSPierre Jolivet {
615e27ab85cSPierre Jolivet   PetscFunctionBegin;
616e27ab85cSPierre Jolivet   C->ops->productsymbolic = MatProductSymbolic_Diagonal_Dense;
617e27ab85cSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
618e27ab85cSPierre Jolivet }
619e27ab85cSPierre Jolivet 
620e27ab85cSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Diagonal_Dense(Mat C)
621e27ab85cSPierre Jolivet {
622e27ab85cSPierre Jolivet   Mat_Product *product = C->product;
623e27ab85cSPierre Jolivet 
624e27ab85cSPierre Jolivet   PetscFunctionBegin;
625e27ab85cSPierre Jolivet   if (product->type == MATPRODUCT_AB || product->type == MATPRODUCT_AtB) PetscCall(MatProductSetFromOptions_Diagonal_Dense_AB(C));
626e27ab85cSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
627e27ab85cSPierre Jolivet }
628e27ab85cSPierre Jolivet 
629345a4b08SToby Isaac /*MC
630345a4b08SToby Isaac    MATDIAGONAL - MATDIAGONAL = "diagonal" - A diagonal matrix type with the diagonal implemented as a `Vec`.  Useful for
631345a4b08SToby Isaac    cases where `VecPointwiseMult()` or `VecPointwiseDivide()` should be thought of as the actions of a linear operator.
632345a4b08SToby Isaac 
633345a4b08SToby Isaac   Level: advanced
634345a4b08SToby Isaac 
635be50c303SSatish Balay .seealso: [](ch_matrices), `Mat`, `MatCreateDiagonal()`, `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()`
636345a4b08SToby Isaac M*/
637345a4b08SToby Isaac PETSC_INTERN PetscErrorCode MatCreate_Diagonal(Mat A)
638345a4b08SToby Isaac {
639345a4b08SToby Isaac   Mat_Diagonal *ctx;
640345a4b08SToby Isaac 
641345a4b08SToby Isaac   PetscFunctionBegin;
642345a4b08SToby Isaac   PetscCall(PetscNew(&ctx));
643345a4b08SToby Isaac   A->data = (void *)ctx;
644345a4b08SToby Isaac 
645345a4b08SToby Isaac   A->structurally_symmetric      = PETSC_BOOL3_TRUE;
646345a4b08SToby Isaac   A->structural_symmetry_eternal = PETSC_TRUE;
647345a4b08SToby Isaac   A->symmetry_eternal            = PETSC_TRUE;
648345a4b08SToby Isaac   A->symmetric                   = PETSC_BOOL3_TRUE;
649345a4b08SToby Isaac   if (!PetscDefined(USE_COMPLEX)) A->hermitian = PETSC_BOOL3_TRUE;
650345a4b08SToby Isaac 
651c3e1b152SPierre Jolivet   A->ops->getrow           = MatGetRow_Diagonal;
652c3e1b152SPierre Jolivet   A->ops->restorerow       = MatRestoreRow_Diagonal;
653345a4b08SToby Isaac   A->ops->mult             = MatMult_Diagonal;
654345a4b08SToby Isaac   A->ops->multadd          = MatMultAdd_Diagonal;
655345a4b08SToby Isaac   A->ops->multtranspose    = MatMult_Diagonal;
656345a4b08SToby Isaac   A->ops->multtransposeadd = MatMultAdd_Diagonal;
657345a4b08SToby Isaac   A->ops->norm             = MatNorm_Diagonal;
658345a4b08SToby Isaac   A->ops->duplicate        = MatDuplicate_Diagonal;
659345a4b08SToby Isaac   A->ops->solve            = MatSolve_Diagonal;
660345a4b08SToby Isaac   A->ops->solvetranspose   = MatSolve_Diagonal;
661345a4b08SToby Isaac   A->ops->shift            = MatShift_Diagonal;
662345a4b08SToby Isaac   A->ops->scale            = MatScale_Diagonal;
663345a4b08SToby Isaac   A->ops->diagonalscale    = MatDiagonalScale_Diagonal;
664345a4b08SToby Isaac   A->ops->getdiagonal      = MatGetDiagonal_Diagonal;
665345a4b08SToby Isaac   A->ops->diagonalset      = MatDiagonalSet_Diagonal;
666345a4b08SToby Isaac   A->ops->view             = MatView_Diagonal;
667345a4b08SToby Isaac   A->ops->zeroentries      = MatZeroEntries_Diagonal;
668345a4b08SToby Isaac   A->ops->destroy          = MatDestroy_Diagonal;
669345a4b08SToby Isaac   A->ops->getinfo          = MatGetInfo_Diagonal;
670345a4b08SToby Isaac   A->ops->axpy             = MatAXPY_Diagonal;
671345a4b08SToby Isaac   A->ops->setup            = MatSetUp_Diagonal;
672*ee8a0a41SPierre Jolivet   A->ops->permute          = MatPermute_Diagonal;
673345a4b08SToby Isaac 
674345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalGetDiagonal_C", MatDiagonalGetDiagonal_Diagonal));
675345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalRestoreDiagonal_C", MatDiagonalRestoreDiagonal_Diagonal));
676345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalGetInverseDiagonal_C", MatDiagonalGetInverseDiagonal_Diagonal));
677345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalRestoreInverseDiagonal_C", MatDiagonalRestoreInverseDiagonal_Diagonal));
678e27ab85cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_diagonal_seqdense_C", MatProductSetFromOptions_Diagonal_Dense));
679e27ab85cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_diagonal_mpidense_C", MatProductSetFromOptions_Diagonal_Dense));
680345a4b08SToby Isaac   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATDIAGONAL));
681345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
682345a4b08SToby Isaac }
683