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; 57*09a7ae42SPierre Jolivet PetscInt rstart = A->rmap->rstart, rend = A->rmap->rend; 58c3e1b152SPierre Jolivet 59c3e1b152SPierre Jolivet PetscFunctionBegin; 60*09a7ae42SPierre Jolivet PetscCheck(row >= rstart && row < rend, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only local rows"); 61c3e1b152SPierre Jolivet if (ncols) *ncols = 1; 62c3e1b152SPierre Jolivet if (cols) { 63c3e1b152SPierre Jolivet if (!mat->col) PetscCall(PetscMalloc1(1, &mat->col)); 64c3e1b152SPierre Jolivet *mat->col = row; 65c3e1b152SPierre Jolivet *cols = mat->col; 66c3e1b152SPierre Jolivet } 67c3e1b152SPierre Jolivet if (vals) { 68c3e1b152SPierre Jolivet const PetscScalar *v; 69c3e1b152SPierre Jolivet 70c3e1b152SPierre Jolivet if (!mat->val) PetscCall(PetscMalloc1(1, &mat->val)); 71c3e1b152SPierre Jolivet PetscCall(VecGetArrayRead(mat->diag, &v)); 72*09a7ae42SPierre Jolivet *mat->val = v[row - rstart]; 73c3e1b152SPierre Jolivet *vals = mat->val; 74c3e1b152SPierre Jolivet PetscCall(VecRestoreArrayRead(mat->diag, &v)); 75c3e1b152SPierre Jolivet } 76c3e1b152SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 77c3e1b152SPierre Jolivet } 78c3e1b152SPierre Jolivet 79345a4b08SToby Isaac static PetscErrorCode MatMult_Diagonal(Mat A, Vec x, Vec y) 80345a4b08SToby Isaac { 81345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 82345a4b08SToby Isaac 83345a4b08SToby Isaac PetscFunctionBegin; 84345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(A)); 85345a4b08SToby Isaac PetscCall(VecPointwiseMult(y, ctx->diag, x)); 86345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 87345a4b08SToby Isaac } 88345a4b08SToby Isaac 89345a4b08SToby Isaac static PetscErrorCode MatMultAdd_Diagonal(Mat mat, Vec v1, Vec v2, Vec v3) 90345a4b08SToby Isaac { 91345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)mat->data; 92345a4b08SToby Isaac 93345a4b08SToby Isaac PetscFunctionBegin; 94345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(mat)); 95345a4b08SToby Isaac if (v2 != v3) { 96345a4b08SToby Isaac PetscCall(VecPointwiseMult(v3, ctx->diag, v1)); 97345a4b08SToby Isaac PetscCall(VecAXPY(v3, 1.0, v2)); 98345a4b08SToby Isaac } else { 99345a4b08SToby Isaac Vec w; 100345a4b08SToby Isaac PetscCall(VecDuplicate(v3, &w)); 101345a4b08SToby Isaac PetscCall(VecPointwiseMult(w, ctx->diag, v1)); 102345a4b08SToby Isaac PetscCall(VecAXPY(v3, 1.0, w)); 103345a4b08SToby Isaac PetscCall(VecDestroy(&w)); 104345a4b08SToby Isaac } 105345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 106345a4b08SToby Isaac } 107345a4b08SToby Isaac 108345a4b08SToby Isaac static PetscErrorCode MatNorm_Diagonal(Mat A, NormType type, PetscReal *nrm) 109345a4b08SToby Isaac { 110345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 111345a4b08SToby Isaac 112345a4b08SToby Isaac PetscFunctionBegin; 113345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(A)); 114345a4b08SToby Isaac type = (type == NORM_FROBENIUS) ? NORM_2 : NORM_INFINITY; 115345a4b08SToby Isaac PetscCall(VecNorm(ctx->diag, type, nrm)); 116345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 117345a4b08SToby Isaac } 118345a4b08SToby Isaac 119345a4b08SToby Isaac static PetscErrorCode MatDuplicate_Diagonal(Mat A, MatDuplicateOption op, Mat *B) 120345a4b08SToby Isaac { 121345a4b08SToby Isaac Mat_Diagonal *actx = (Mat_Diagonal *)A->data; 122345a4b08SToby Isaac 123345a4b08SToby Isaac PetscFunctionBegin; 124345a4b08SToby Isaac PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B)); 125345a4b08SToby Isaac PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 126345a4b08SToby Isaac PetscCall(MatSetBlockSizesFromMats(*B, A, A)); 127345a4b08SToby Isaac PetscCall(MatSetType(*B, MATDIAGONAL)); 128345a4b08SToby Isaac PetscCall(PetscLayoutReference(A->rmap, &(*B)->rmap)); 129345a4b08SToby Isaac PetscCall(PetscLayoutReference(A->cmap, &(*B)->cmap)); 130345a4b08SToby Isaac if (op == MAT_COPY_VALUES) { 131345a4b08SToby Isaac Mat_Diagonal *bctx = (Mat_Diagonal *)(*B)->data; 132345a4b08SToby Isaac 133345a4b08SToby Isaac PetscCall(MatSetUp(A)); 134345a4b08SToby Isaac PetscCall(MatSetUp(*B)); 135345a4b08SToby Isaac PetscCall(VecCopy(actx->diag, bctx->diag)); 136345a4b08SToby Isaac PetscCall(VecCopy(actx->inv_diag, bctx->inv_diag)); 137345a4b08SToby Isaac bctx->diag_valid = actx->diag_valid; 138345a4b08SToby Isaac bctx->inv_diag_valid = actx->inv_diag_valid; 139345a4b08SToby Isaac } 140345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 141345a4b08SToby Isaac } 142345a4b08SToby Isaac 143345a4b08SToby Isaac /*@ 144345a4b08SToby Isaac MatDiagonalGetDiagonal - Get the diagonal of a `MATDIAGONAL` 145345a4b08SToby Isaac 146345a4b08SToby Isaac Input Parameter: 147345a4b08SToby Isaac . A - the `MATDIAGONAL` 148345a4b08SToby Isaac 149345a4b08SToby Isaac Output Parameter: 150345a4b08SToby Isaac . diag - the `Vec` that defines the diagonal 151345a4b08SToby Isaac 152345a4b08SToby Isaac Level: developer 153345a4b08SToby Isaac 154345a4b08SToby Isaac Note: 155345a4b08SToby Isaac The user must call 156345a4b08SToby Isaac `MatDiagonalRestoreDiagonal()` before using the matrix again. 157345a4b08SToby Isaac 158345a4b08SToby Isaac For a copy of the diagonal values, rather than a reference, use `MatGetDiagonal()` 159345a4b08SToby Isaac 160345a4b08SToby 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()` 161345a4b08SToby Isaac 162be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()`, `MatGetDiagonal()` 163345a4b08SToby Isaac @*/ 164345a4b08SToby Isaac PetscErrorCode MatDiagonalGetDiagonal(Mat A, Vec *diag) 165345a4b08SToby Isaac { 166345a4b08SToby Isaac PetscFunctionBegin; 167345a4b08SToby Isaac PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 1684f572ea9SToby Isaac PetscAssertPointer(diag, 2); 169345a4b08SToby Isaac *diag = NULL; 170835f2295SStefano Zampini PetscUseMethod(A, "MatDiagonalGetDiagonal_C", (Mat, Vec *), (A, diag)); 171345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 172345a4b08SToby Isaac } 173345a4b08SToby Isaac 174345a4b08SToby Isaac static PetscErrorCode MatDiagonalGetDiagonal_Diagonal(Mat A, Vec *diag) 175345a4b08SToby Isaac { 176345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 177345a4b08SToby Isaac 178345a4b08SToby Isaac PetscFunctionBegin; 179345a4b08SToby Isaac PetscCall(MatSetUp(A)); 180345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(A)); 181345a4b08SToby Isaac *diag = ctx->diag; 182345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)*diag, &ctx->diag_state)); 183345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 184345a4b08SToby Isaac } 185345a4b08SToby Isaac 186345a4b08SToby Isaac /*@ 187345a4b08SToby Isaac MatDiagonalRestoreDiagonal - Restore the diagonal of a `MATDIAGONAL` 188345a4b08SToby Isaac 189345a4b08SToby Isaac Input Parameters: 190345a4b08SToby Isaac + A - the `MATDIAGONAL` 191345a4b08SToby Isaac - diag - the `Vec` obtained from `MatDiagonalGetDiagonal()` 192345a4b08SToby Isaac 193345a4b08SToby Isaac Level: developer 194345a4b08SToby Isaac 195345a4b08SToby Isaac Note: 196345a4b08SToby Isaac Use `MatDiagonalSet()` to change the values by copy, rather than reference. 197345a4b08SToby Isaac 198be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalGetDiagonal()` 199345a4b08SToby Isaac @*/ 200345a4b08SToby Isaac PetscErrorCode MatDiagonalRestoreDiagonal(Mat A, Vec *diag) 201345a4b08SToby Isaac { 202345a4b08SToby Isaac PetscFunctionBegin; 203345a4b08SToby Isaac PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2044f572ea9SToby Isaac PetscAssertPointer(diag, 2); 205835f2295SStefano Zampini PetscUseMethod(A, "MatDiagonalRestoreDiagonal_C", (Mat, Vec *), (A, diag)); 206345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 207345a4b08SToby Isaac } 208345a4b08SToby Isaac 209345a4b08SToby Isaac static PetscErrorCode MatDiagonalRestoreDiagonal_Diagonal(Mat A, Vec *diag) 210345a4b08SToby Isaac { 211345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 212345a4b08SToby Isaac PetscObjectState diag_state; 213345a4b08SToby Isaac 214345a4b08SToby Isaac PetscFunctionBegin; 215345a4b08SToby Isaac PetscCheck(ctx->diag == *diag, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Restored a different diagonal vector"); 216345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 217345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)*diag, &diag_state)); 218345a4b08SToby Isaac if (ctx->diag_state != diag_state) { 219345a4b08SToby Isaac PetscCall(PetscObjectStateIncrease((PetscObject)A)); 220345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 221345a4b08SToby Isaac } 222345a4b08SToby Isaac *diag = NULL; 223345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 224345a4b08SToby Isaac } 225345a4b08SToby Isaac 226345a4b08SToby Isaac /*@ 227345a4b08SToby Isaac MatDiagonalGetInverseDiagonal - Get the inverse diagonal of a `MATDIAGONAL` 228345a4b08SToby Isaac 229345a4b08SToby Isaac Input Parameter: 230345a4b08SToby Isaac . A - the `MATDIAGONAL` 231345a4b08SToby Isaac 232345a4b08SToby Isaac Output Parameter: 233345a4b08SToby Isaac . inv_diag - the `Vec` that defines the inverse diagonal 234345a4b08SToby Isaac 235345a4b08SToby Isaac Level: developer 236345a4b08SToby Isaac 237345a4b08SToby Isaac Note: 238345a4b08SToby Isaac The user must call 239345a4b08SToby Isaac `MatDiagonalRestoreInverseDiagonal()` before using the matrix again. 240345a4b08SToby Isaac 241345a4b08SToby Isaac If a matrix is created only to call `MatSolve()` (which happens for `MATLMVMDIAGBROYDEN`), 242345a4b08SToby Isaac using `MatDiagonalGetInverseDiagonal()` and `MatDiagonalRestoreInverseDiagonal()` avoids copies 243345a4b08SToby Isaac and avoids any call to `VecReciprocal()`. 244345a4b08SToby Isaac 245be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MATLMVMBROYDEN`, `MatSolve()` 246345a4b08SToby Isaac @*/ 247345a4b08SToby Isaac PetscErrorCode MatDiagonalGetInverseDiagonal(Mat A, Vec *inv_diag) 248345a4b08SToby Isaac { 249345a4b08SToby Isaac PetscFunctionBegin; 250345a4b08SToby Isaac PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2514f572ea9SToby Isaac PetscAssertPointer(inv_diag, 2); 252345a4b08SToby Isaac *inv_diag = NULL; 253835f2295SStefano Zampini PetscUseMethod(A, "MatDiagonalGetInverseDiagonal_C", (Mat, Vec *), (A, inv_diag)); 254345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 255345a4b08SToby Isaac } 256345a4b08SToby Isaac 257345a4b08SToby Isaac static PetscErrorCode MatDiagonalGetInverseDiagonal_Diagonal(Mat A, Vec *inv_diag) 258345a4b08SToby Isaac { 259345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 260345a4b08SToby Isaac 261345a4b08SToby Isaac PetscFunctionBegin; 262345a4b08SToby Isaac PetscCall(MatSetUp(A)); 263345a4b08SToby Isaac PetscCall(MatDiagonalSetUpInverseDiagonal(A)); 264345a4b08SToby Isaac *inv_diag = ctx->inv_diag; 265345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)*inv_diag, &ctx->inv_diag_state)); 266345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 267345a4b08SToby Isaac } 268345a4b08SToby Isaac 269345a4b08SToby Isaac /*@ 270345a4b08SToby Isaac MatDiagonalRestoreInverseDiagonal - Restore the inverse diagonal of a `MATDIAGONAL` 271345a4b08SToby Isaac 272345a4b08SToby Isaac Input Parameters: 273345a4b08SToby Isaac + A - the `MATDIAGONAL` 274345a4b08SToby Isaac - inv_diag - the `Vec` obtained from `MatDiagonalGetInverseDiagonal()` 275345a4b08SToby Isaac 276345a4b08SToby Isaac Level: developer 277345a4b08SToby Isaac 278be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalGetInverseDiagonal()` 279345a4b08SToby Isaac @*/ 280345a4b08SToby Isaac PetscErrorCode MatDiagonalRestoreInverseDiagonal(Mat A, Vec *inv_diag) 281345a4b08SToby Isaac { 282345a4b08SToby Isaac PetscFunctionBegin; 283345a4b08SToby Isaac PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2844f572ea9SToby Isaac PetscAssertPointer(inv_diag, 2); 285835f2295SStefano Zampini PetscUseMethod(A, "MatDiagonalRestoreInverseDiagonal_C", (Mat, Vec *), (A, inv_diag)); 286345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 287345a4b08SToby Isaac } 288345a4b08SToby Isaac 289345a4b08SToby Isaac static PetscErrorCode MatDiagonalRestoreInverseDiagonal_Diagonal(Mat A, Vec *inv_diag) 290345a4b08SToby Isaac { 291345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 292345a4b08SToby Isaac PetscObjectState inv_diag_state; 293345a4b08SToby Isaac 294345a4b08SToby Isaac PetscFunctionBegin; 295345a4b08SToby Isaac PetscCheck(ctx->inv_diag == *inv_diag, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Restored a different diagonal vector"); 296345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_TRUE; 297345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)*inv_diag, &inv_diag_state)); 298345a4b08SToby Isaac if (ctx->inv_diag_state != inv_diag_state) { 299345a4b08SToby Isaac PetscCall(PetscObjectStateIncrease((PetscObject)A)); 300345a4b08SToby Isaac ctx->diag_valid = PETSC_FALSE; 301345a4b08SToby Isaac } 302345a4b08SToby Isaac *inv_diag = NULL; 303345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 304345a4b08SToby Isaac } 305345a4b08SToby Isaac 306ee8a0a41SPierre Jolivet static PetscErrorCode MatPermute_Diagonal(Mat A, IS rowp, IS colp, Mat *B) 307ee8a0a41SPierre Jolivet { 308ee8a0a41SPierre Jolivet Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 309ee8a0a41SPierre Jolivet Vec v; 310ee8a0a41SPierre Jolivet 311ee8a0a41SPierre Jolivet PetscFunctionBegin; 312ee8a0a41SPierre Jolivet PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same"); 313ee8a0a41SPierre Jolivet PetscCall(VecDuplicate(ctx->diag, &v)); 314ee8a0a41SPierre Jolivet PetscCall(VecCopy(ctx->diag, v)); 315ee8a0a41SPierre Jolivet PetscCall(VecPermute(v, rowp, PETSC_FALSE)); 316ee8a0a41SPierre Jolivet PetscCall(MatCreateDiagonal(v, B)); 317ee8a0a41SPierre Jolivet PetscCall(VecDestroy(&v)); 318ee8a0a41SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 319ee8a0a41SPierre Jolivet } 320ee8a0a41SPierre Jolivet 3218a1df862SToby Isaac static PetscErrorCode MatSetRandom_Diagonal(Mat A, PetscRandom rand) 3228a1df862SToby Isaac { 3238a1df862SToby Isaac Vec d; 3248a1df862SToby Isaac 3258a1df862SToby Isaac PetscFunctionBegin; 3268a1df862SToby Isaac PetscCall(MatDiagonalGetDiagonal(A, &d)); 3278a1df862SToby Isaac PetscCall(VecSetRandom(d, rand)); 3288a1df862SToby Isaac PetscCall(MatDiagonalRestoreDiagonal(A, &d)); 3298a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 3308a1df862SToby Isaac } 3318a1df862SToby Isaac 332345a4b08SToby Isaac static PetscErrorCode MatDestroy_Diagonal(Mat mat) 333345a4b08SToby Isaac { 334345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)mat->data; 335345a4b08SToby Isaac 336345a4b08SToby Isaac PetscFunctionBegin; 337345a4b08SToby Isaac PetscCall(VecDestroy(&ctx->diag)); 338345a4b08SToby Isaac PetscCall(VecDestroy(&ctx->inv_diag)); 339c3e1b152SPierre Jolivet PetscCall(PetscFree(ctx->col)); 340c3e1b152SPierre Jolivet PetscCall(PetscFree(ctx->val)); 341345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalGetDiagonal_C", NULL)); 342345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalRestoreDiagonal_C", NULL)); 343345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalGetInverseDiagonal_C", NULL)); 344345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalRestoreInverseDiagonal_C", NULL)); 345e27ab85cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_diagonal_seqdense_C", NULL)); 346e27ab85cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatProductSetFromOptions_diagonal_mpidense_C", NULL)); 347345a4b08SToby Isaac PetscCall(PetscFree(mat->data)); 348345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 349345a4b08SToby Isaac } 350345a4b08SToby Isaac 351345a4b08SToby Isaac static PetscErrorCode MatView_Diagonal(Mat J, PetscViewer viewer) 352345a4b08SToby Isaac { 353345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)J->data; 3549f196a02SMartin Diehl PetscBool isascii; 355345a4b08SToby Isaac 356345a4b08SToby Isaac PetscFunctionBegin; 3579f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 3589f196a02SMartin Diehl if (isascii) { 359345a4b08SToby Isaac PetscViewerFormat format; 360345a4b08SToby Isaac 361345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 3628a1df862SToby Isaac PetscCall(PetscViewerGetFormat(viewer, &format)); 3638a1df862SToby Isaac if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 3648a1df862SToby Isaac PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ctx->diag, viewer)); 3658a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 3668a1df862SToby Isaac } 367345a4b08SToby Isaac PetscCall(VecView(ctx->diag, viewer)); 368345a4b08SToby Isaac } 369345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 370345a4b08SToby Isaac } 371345a4b08SToby Isaac 372345a4b08SToby Isaac static PetscErrorCode MatGetDiagonal_Diagonal(Mat J, Vec x) 373345a4b08SToby Isaac { 374345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)J->data; 375345a4b08SToby Isaac 376345a4b08SToby Isaac PetscFunctionBegin; 377345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 378345a4b08SToby Isaac PetscCall(VecCopy(ctx->diag, x)); 379345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 380345a4b08SToby Isaac } 381345a4b08SToby Isaac 382345a4b08SToby Isaac static PetscErrorCode MatDiagonalSet_Diagonal(Mat J, Vec D, InsertMode is) 383345a4b08SToby Isaac { 384345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)J->data; 385345a4b08SToby Isaac 386345a4b08SToby Isaac PetscFunctionBegin; 387345a4b08SToby Isaac switch (is) { 388345a4b08SToby Isaac case ADD_VALUES: 389345a4b08SToby Isaac case ADD_ALL_VALUES: 390345a4b08SToby Isaac case ADD_BC_VALUES: 391345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 392345a4b08SToby Isaac PetscCall(VecAXPY(ctx->diag, 1.0, D)); 393345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 394345a4b08SToby Isaac break; 395345a4b08SToby Isaac case INSERT_VALUES: 396345a4b08SToby Isaac case INSERT_BC_VALUES: 397345a4b08SToby Isaac case INSERT_ALL_VALUES: 398345a4b08SToby Isaac PetscCall(MatSetUp(J)); 399345a4b08SToby Isaac PetscCall(VecCopy(D, ctx->diag)); 400345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 401345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 402345a4b08SToby Isaac break; 403345a4b08SToby Isaac case MAX_VALUES: 404345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 405345a4b08SToby Isaac PetscCall(VecPointwiseMax(ctx->diag, D, ctx->diag)); 406345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 407345a4b08SToby Isaac break; 408345a4b08SToby Isaac case MIN_VALUES: 409345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 410345a4b08SToby Isaac PetscCall(VecPointwiseMin(ctx->diag, D, ctx->diag)); 411345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 412345a4b08SToby Isaac break; 413345a4b08SToby Isaac case NOT_SET_VALUES: 414345a4b08SToby Isaac break; 415345a4b08SToby Isaac } 416345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 417345a4b08SToby Isaac } 418345a4b08SToby Isaac 419345a4b08SToby Isaac static PetscErrorCode MatShift_Diagonal(Mat Y, PetscScalar a) 420345a4b08SToby Isaac { 421345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data; 422345a4b08SToby Isaac 423345a4b08SToby Isaac PetscFunctionBegin; 424345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(Y)); 425345a4b08SToby Isaac PetscCall(VecShift(ctx->diag, a)); 426345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 427345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 428345a4b08SToby Isaac } 429345a4b08SToby Isaac 430345a4b08SToby Isaac static PetscErrorCode MatScale_Diagonal(Mat Y, PetscScalar a) 431345a4b08SToby Isaac { 432345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data; 433345a4b08SToby Isaac 434345a4b08SToby Isaac PetscFunctionBegin; 435345a4b08SToby Isaac PetscCall(MatSetUp(Y)); 436345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(Y)); 437345a4b08SToby Isaac PetscCall(VecScale(ctx->diag, a)); 438345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 439345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 440345a4b08SToby Isaac } 441345a4b08SToby Isaac 442345a4b08SToby Isaac static PetscErrorCode MatDiagonalScale_Diagonal(Mat Y, Vec l, Vec r) 443345a4b08SToby Isaac { 444345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data; 445345a4b08SToby Isaac 446345a4b08SToby Isaac PetscFunctionBegin; 447345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(Y)); 448345a4b08SToby Isaac if (l) { 449345a4b08SToby Isaac PetscCall(VecPointwiseMult(ctx->diag, ctx->diag, l)); 450345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 451345a4b08SToby Isaac } 452345a4b08SToby Isaac if (r) { 453345a4b08SToby Isaac PetscCall(VecPointwiseMult(ctx->diag, ctx->diag, r)); 454345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 455345a4b08SToby Isaac } 456345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 457345a4b08SToby Isaac } 458345a4b08SToby Isaac 4598a1df862SToby Isaac static PetscErrorCode MatConjugate_Diagonal(Mat Y) 4608a1df862SToby Isaac { 4618a1df862SToby Isaac PetscFunctionBegin; 4628a1df862SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data; 4638a1df862SToby Isaac 4648a1df862SToby Isaac if (ctx->diag_valid) { 4658a1df862SToby Isaac PetscCall(VecConjugate(ctx->diag)); 4668a1df862SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)ctx->diag, &ctx->diag_state)); 4678a1df862SToby Isaac } 4688a1df862SToby Isaac if (ctx->inv_diag_valid) { 4698a1df862SToby Isaac PetscCall(VecConjugate(ctx->inv_diag)); 4708a1df862SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)ctx->inv_diag, &ctx->inv_diag_state)); 4718a1df862SToby Isaac } 4728a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 4738a1df862SToby Isaac } 4748a1df862SToby Isaac 4758a1df862SToby Isaac static PetscErrorCode MatTranspose_Diagonal(Mat A, MatReuse reuse, Mat *matout) 4768a1df862SToby Isaac { 4778a1df862SToby Isaac PetscFunctionBegin; 4788a1df862SToby Isaac if (reuse == MAT_INPLACE_MATRIX) { 4798a1df862SToby Isaac PetscLayout tmplayout = A->rmap; 4808a1df862SToby Isaac 4818a1df862SToby Isaac A->rmap = A->cmap; 4828a1df862SToby Isaac A->cmap = tmplayout; 4838a1df862SToby Isaac } else { 4848a1df862SToby Isaac Vec diag, newdiag; 4858a1df862SToby Isaac if (reuse == MAT_INITIAL_MATRIX) { 4868a1df862SToby Isaac PetscCall(MatDiagonalGetDiagonal(A, &diag)); 4878a1df862SToby Isaac PetscCall(VecDuplicate(diag, &newdiag)); 4888a1df862SToby Isaac PetscCall(VecCopy(diag, newdiag)); 4898a1df862SToby Isaac PetscCall(MatDiagonalRestoreDiagonal(A, &diag)); 4908a1df862SToby Isaac PetscCall(MatCreateDiagonal(newdiag, matout)); 4918a1df862SToby Isaac PetscCall(VecDestroy(&newdiag)); 4928a1df862SToby Isaac } else { 4938a1df862SToby Isaac PetscCall(MatDiagonalGetDiagonal(A, &diag)); 4948a1df862SToby Isaac PetscCall(MatDiagonalGetDiagonal(*matout, &newdiag)); 4958a1df862SToby Isaac PetscCall(VecCopy(diag, newdiag)); 4968a1df862SToby Isaac PetscCall(MatDiagonalRestoreDiagonal(*matout, &newdiag)); 4978a1df862SToby Isaac PetscCall(MatDiagonalRestoreDiagonal(A, &diag)); 4988a1df862SToby Isaac } 4998a1df862SToby Isaac } 5008a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 5018a1df862SToby Isaac } 5028a1df862SToby Isaac 503345a4b08SToby Isaac static PetscErrorCode MatSetUp_Diagonal(Mat A) 504345a4b08SToby Isaac { 505345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 506345a4b08SToby Isaac 507345a4b08SToby Isaac PetscFunctionBegin; 508345a4b08SToby Isaac if (!ctx->diag) { 509345a4b08SToby Isaac PetscCall(PetscLayoutSetUp(A->cmap)); 510345a4b08SToby Isaac PetscCall(PetscLayoutSetUp(A->rmap)); 511345a4b08SToby Isaac PetscCall(MatCreateVecs(A, &ctx->diag, NULL)); 512345a4b08SToby Isaac PetscCall(VecDuplicate(ctx->diag, &ctx->inv_diag)); 513345a4b08SToby Isaac PetscCall(VecZeroEntries(ctx->diag)); 514345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 515345a4b08SToby Isaac } 516345a4b08SToby Isaac A->assembled = PETSC_TRUE; 517345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 518345a4b08SToby Isaac } 519345a4b08SToby Isaac 520345a4b08SToby Isaac static PetscErrorCode MatZeroEntries_Diagonal(Mat Y) 521345a4b08SToby Isaac { 522345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data; 523345a4b08SToby Isaac 524345a4b08SToby Isaac PetscFunctionBegin; 525345a4b08SToby Isaac PetscCall(MatSetUp(Y)); 526345a4b08SToby Isaac PetscCall(VecZeroEntries(ctx->diag)); 527345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 528345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 529345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 530345a4b08SToby Isaac } 531345a4b08SToby Isaac 53266976f2fSJacob Faibussowitsch static PetscErrorCode MatSolve_Diagonal(Mat matin, Vec b, Vec x) 533345a4b08SToby Isaac { 534345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)matin->data; 535345a4b08SToby Isaac 536345a4b08SToby Isaac PetscFunctionBegin; 537345a4b08SToby Isaac PetscCall(MatDiagonalSetUpInverseDiagonal(matin)); 538345a4b08SToby Isaac PetscCall(VecPointwiseMult(x, b, ctx->inv_diag)); 539345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 540345a4b08SToby Isaac } 541345a4b08SToby Isaac 54266976f2fSJacob Faibussowitsch static PetscErrorCode MatGetInfo_Diagonal(Mat A, MatInfoType flag, MatInfo *info) 543345a4b08SToby Isaac { 544345a4b08SToby Isaac PetscFunctionBegin; 545345a4b08SToby Isaac info->block_size = 1.0; 546345a4b08SToby Isaac info->nz_allocated = A->cmap->N; 547345a4b08SToby Isaac info->nz_used = A->cmap->N; 548345a4b08SToby Isaac info->nz_unneeded = 0.0; 549345a4b08SToby Isaac info->assemblies = A->num_ass; 550345a4b08SToby Isaac info->mallocs = 0.0; 551345a4b08SToby Isaac info->memory = 0; 552345a4b08SToby Isaac info->fill_ratio_given = 0; 553345a4b08SToby Isaac info->fill_ratio_needed = 0; 554345a4b08SToby Isaac info->factor_mallocs = 0; 555345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 556345a4b08SToby Isaac } 557345a4b08SToby Isaac 558345a4b08SToby Isaac /*@ 559345a4b08SToby Isaac MatCreateDiagonal - Creates a matrix defined by a given vector along its diagonal. 560345a4b08SToby Isaac 561345a4b08SToby Isaac Collective 562345a4b08SToby Isaac 563345a4b08SToby Isaac Input Parameter: 564345a4b08SToby Isaac . diag - vector for the diagonal 565345a4b08SToby Isaac 566345a4b08SToby Isaac Output Parameter: 567345a4b08SToby Isaac . J - the diagonal matrix 568345a4b08SToby Isaac 569345a4b08SToby Isaac Level: advanced 570345a4b08SToby Isaac 571345a4b08SToby Isaac Notes: 572345a4b08SToby Isaac Only supports square matrices with the same number of local rows and columns. 573345a4b08SToby Isaac 574345a4b08SToby Isaac The input vector `diag` will be referenced internally: any changes to `diag` 575345a4b08SToby Isaac will affect the matrix `J`. 576345a4b08SToby Isaac 577be50c303SSatish Balay .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatSolve()` 578345a4b08SToby Isaac `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()` 579345a4b08SToby Isaac @*/ 580345a4b08SToby Isaac PetscErrorCode MatCreateDiagonal(Vec diag, Mat *J) 581345a4b08SToby Isaac { 582345a4b08SToby Isaac PetscFunctionBegin; 583345a4b08SToby Isaac PetscValidHeaderSpecific(diag, VEC_CLASSID, 1); 584345a4b08SToby Isaac PetscCall(MatCreate(PetscObjectComm((PetscObject)diag), J)); 585345a4b08SToby Isaac PetscInt m, M; 586345a4b08SToby Isaac PetscCall(VecGetLocalSize(diag, &m)); 587345a4b08SToby Isaac PetscCall(VecGetSize(diag, &M)); 588345a4b08SToby Isaac PetscCall(MatSetSizes(*J, m, m, M, M)); 589345a4b08SToby Isaac PetscCall(MatSetType(*J, MATDIAGONAL)); 590345a4b08SToby Isaac 591345a4b08SToby Isaac PetscLayout map; 592345a4b08SToby Isaac PetscCall(VecGetLayout(diag, &map)); 593345a4b08SToby Isaac PetscCall(MatSetLayouts(*J, map, map)); 594345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)(*J)->data; 595345a4b08SToby Isaac PetscCall(PetscObjectReference((PetscObject)diag)); 596345a4b08SToby Isaac PetscCall(VecDestroy(&ctx->diag)); 597345a4b08SToby Isaac PetscCall(VecDestroy(&ctx->inv_diag)); 598345a4b08SToby Isaac ctx->diag = diag; 599345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 600345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 601345a4b08SToby Isaac VecType type; 602345a4b08SToby Isaac PetscCall(VecDuplicate(ctx->diag, &ctx->inv_diag)); 603345a4b08SToby Isaac PetscCall(VecGetType(diag, &type)); 604345a4b08SToby Isaac PetscCall(PetscFree((*J)->defaultvectype)); 605345a4b08SToby Isaac PetscCall(PetscStrallocpy(type, &(*J)->defaultvectype)); 606345a4b08SToby Isaac PetscCall(MatSetUp(*J)); 607c3e1b152SPierre Jolivet ctx->col = NULL; 608c3e1b152SPierre Jolivet ctx->val = NULL; 609345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 610345a4b08SToby Isaac } 611345a4b08SToby Isaac 612e27ab85cSPierre Jolivet static PetscErrorCode MatProductNumeric_Diagonal_Dense(Mat C) 613e27ab85cSPierre Jolivet { 614e27ab85cSPierre Jolivet Mat A, B; 615e27ab85cSPierre Jolivet Mat_Diagonal *a; 616e27ab85cSPierre Jolivet PetscScalar *c; 617e27ab85cSPierre Jolivet const PetscScalar *b, *alpha; 618e27ab85cSPierre Jolivet PetscInt ldb, ldc; 619e27ab85cSPierre Jolivet 620e27ab85cSPierre Jolivet PetscFunctionBegin; 621e27ab85cSPierre Jolivet MatCheckProduct(C, 1); 622e27ab85cSPierre Jolivet A = C->product->A; 623e27ab85cSPierre Jolivet B = C->product->B; 624e27ab85cSPierre Jolivet a = (Mat_Diagonal *)A->data; 625e27ab85cSPierre Jolivet PetscCall(VecGetArrayRead(a->diag, &alpha)); 626e27ab85cSPierre Jolivet PetscCall(MatDenseGetLDA(B, &ldb)); 627e27ab85cSPierre Jolivet PetscCall(MatDenseGetLDA(C, &ldc)); 628e27ab85cSPierre Jolivet PetscCall(MatDenseGetArrayRead(B, &b)); 629e27ab85cSPierre Jolivet PetscCall(MatDenseGetArrayWrite(C, &c)); 630e27ab85cSPierre Jolivet for (PetscInt j = 0; j < B->cmap->N; j++) 631e27ab85cSPierre Jolivet for (PetscInt i = 0; i < B->rmap->n; i++) c[i + j * ldc] = alpha[i] * b[i + j * ldb]; 632e27ab85cSPierre Jolivet PetscCall(MatDenseRestoreArrayWrite(C, &c)); 633e27ab85cSPierre Jolivet PetscCall(MatDenseRestoreArrayRead(B, &b)); 634e27ab85cSPierre Jolivet PetscCall(VecRestoreArrayRead(a->diag, &alpha)); 635e27ab85cSPierre Jolivet PetscCall(PetscLogFlops(B->cmap->N * B->rmap->n)); 636e27ab85cSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 637e27ab85cSPierre Jolivet } 638e27ab85cSPierre Jolivet 639e27ab85cSPierre Jolivet static PetscErrorCode MatProductSymbolic_Diagonal_Dense(Mat C) 640e27ab85cSPierre Jolivet { 641e27ab85cSPierre Jolivet Mat A, B; 642e27ab85cSPierre Jolivet PetscInt n, N, m, M; 643e27ab85cSPierre Jolivet 644e27ab85cSPierre Jolivet PetscFunctionBegin; 645e27ab85cSPierre Jolivet MatCheckProduct(C, 1); 646e27ab85cSPierre Jolivet PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty"); 647e27ab85cSPierre Jolivet A = C->product->A; 648e27ab85cSPierre Jolivet B = C->product->B; 649e27ab85cSPierre Jolivet PetscCall(MatDiagonalSetUpDiagonal(A)); 650e27ab85cSPierre Jolivet PetscCall(MatGetLocalSize(C, &m, &n)); 651e27ab85cSPierre Jolivet PetscCall(MatGetSize(C, &M, &N)); 652e27ab85cSPierre Jolivet if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) { 653e27ab85cSPierre Jolivet PetscCall(MatGetLocalSize(B, NULL, &n)); 654e27ab85cSPierre Jolivet PetscCall(MatGetSize(B, NULL, &N)); 655e27ab85cSPierre Jolivet PetscCall(MatGetLocalSize(A, &m, NULL)); 656e27ab85cSPierre Jolivet PetscCall(MatGetSize(A, &M, NULL)); 657e27ab85cSPierre Jolivet PetscCall(MatSetSizes(C, m, n, M, N)); 658e27ab85cSPierre Jolivet } 659e27ab85cSPierre Jolivet PetscCall(MatSetType(C, ((PetscObject)B)->type_name)); 660e27ab85cSPierre Jolivet PetscCall(MatSetUp(C)); 661e27ab85cSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Diagonal_Dense; 662e27ab85cSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 663e27ab85cSPierre Jolivet } 664e27ab85cSPierre Jolivet 665e27ab85cSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Diagonal_Dense_AB(Mat C) 666e27ab85cSPierre Jolivet { 667e27ab85cSPierre Jolivet PetscFunctionBegin; 668e27ab85cSPierre Jolivet C->ops->productsymbolic = MatProductSymbolic_Diagonal_Dense; 669e27ab85cSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 670e27ab85cSPierre Jolivet } 671e27ab85cSPierre Jolivet 672e27ab85cSPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Diagonal_Dense(Mat C) 673e27ab85cSPierre Jolivet { 674e27ab85cSPierre Jolivet Mat_Product *product = C->product; 675e27ab85cSPierre Jolivet 676e27ab85cSPierre Jolivet PetscFunctionBegin; 677e27ab85cSPierre Jolivet if (product->type == MATPRODUCT_AB || product->type == MATPRODUCT_AtB) PetscCall(MatProductSetFromOptions_Diagonal_Dense_AB(C)); 678e27ab85cSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 679e27ab85cSPierre Jolivet } 680e27ab85cSPierre Jolivet 681345a4b08SToby Isaac /*MC 682345a4b08SToby Isaac MATDIAGONAL - MATDIAGONAL = "diagonal" - A diagonal matrix type with the diagonal implemented as a `Vec`. Useful for 683345a4b08SToby Isaac cases where `VecPointwiseMult()` or `VecPointwiseDivide()` should be thought of as the actions of a linear operator. 684345a4b08SToby Isaac 685345a4b08SToby Isaac Level: advanced 686345a4b08SToby Isaac 687be50c303SSatish Balay .seealso: [](ch_matrices), `Mat`, `MatCreateDiagonal()`, `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()` 688345a4b08SToby Isaac M*/ 689345a4b08SToby Isaac PETSC_INTERN PetscErrorCode MatCreate_Diagonal(Mat A) 690345a4b08SToby Isaac { 691345a4b08SToby Isaac Mat_Diagonal *ctx; 692345a4b08SToby Isaac 693345a4b08SToby Isaac PetscFunctionBegin; 694345a4b08SToby Isaac PetscCall(PetscNew(&ctx)); 695345a4b08SToby Isaac A->data = (void *)ctx; 696345a4b08SToby Isaac 697345a4b08SToby Isaac A->structurally_symmetric = PETSC_BOOL3_TRUE; 698345a4b08SToby Isaac A->structural_symmetry_eternal = PETSC_TRUE; 699345a4b08SToby Isaac A->symmetry_eternal = PETSC_TRUE; 700345a4b08SToby Isaac A->symmetric = PETSC_BOOL3_TRUE; 701345a4b08SToby Isaac if (!PetscDefined(USE_COMPLEX)) A->hermitian = PETSC_BOOL3_TRUE; 702345a4b08SToby Isaac 703c3e1b152SPierre Jolivet A->ops->getrow = MatGetRow_Diagonal; 704345a4b08SToby Isaac A->ops->mult = MatMult_Diagonal; 705345a4b08SToby Isaac A->ops->multadd = MatMultAdd_Diagonal; 706345a4b08SToby Isaac A->ops->multtranspose = MatMult_Diagonal; 707345a4b08SToby Isaac A->ops->multtransposeadd = MatMultAdd_Diagonal; 708345a4b08SToby Isaac A->ops->norm = MatNorm_Diagonal; 709345a4b08SToby Isaac A->ops->duplicate = MatDuplicate_Diagonal; 710345a4b08SToby Isaac A->ops->solve = MatSolve_Diagonal; 711345a4b08SToby Isaac A->ops->solvetranspose = MatSolve_Diagonal; 712345a4b08SToby Isaac A->ops->shift = MatShift_Diagonal; 713345a4b08SToby Isaac A->ops->scale = MatScale_Diagonal; 714345a4b08SToby Isaac A->ops->diagonalscale = MatDiagonalScale_Diagonal; 715345a4b08SToby Isaac A->ops->getdiagonal = MatGetDiagonal_Diagonal; 716345a4b08SToby Isaac A->ops->diagonalset = MatDiagonalSet_Diagonal; 717345a4b08SToby Isaac A->ops->view = MatView_Diagonal; 718345a4b08SToby Isaac A->ops->zeroentries = MatZeroEntries_Diagonal; 719345a4b08SToby Isaac A->ops->destroy = MatDestroy_Diagonal; 720345a4b08SToby Isaac A->ops->getinfo = MatGetInfo_Diagonal; 721345a4b08SToby Isaac A->ops->axpy = MatAXPY_Diagonal; 722345a4b08SToby Isaac A->ops->setup = MatSetUp_Diagonal; 723ee8a0a41SPierre Jolivet A->ops->permute = MatPermute_Diagonal; 7248a1df862SToby Isaac A->ops->setrandom = MatSetRandom_Diagonal; 7258a1df862SToby Isaac A->ops->conjugate = MatConjugate_Diagonal; 7268a1df862SToby Isaac A->ops->transpose = MatTranspose_Diagonal; 727345a4b08SToby Isaac 728345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalGetDiagonal_C", MatDiagonalGetDiagonal_Diagonal)); 729345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalRestoreDiagonal_C", MatDiagonalRestoreDiagonal_Diagonal)); 730345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalGetInverseDiagonal_C", MatDiagonalGetInverseDiagonal_Diagonal)); 731345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalRestoreInverseDiagonal_C", MatDiagonalRestoreInverseDiagonal_Diagonal)); 732e27ab85cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_diagonal_seqdense_C", MatProductSetFromOptions_Diagonal_Dense)); 733e27ab85cSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_diagonal_mpidense_C", MatProductSetFromOptions_Diagonal_Dense)); 734345a4b08SToby Isaac PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATDIAGONAL)); 735345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 736345a4b08SToby Isaac } 737