1345a4b08SToby Isaac 2345a4b08SToby Isaac #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3345a4b08SToby Isaac 4345a4b08SToby Isaac typedef struct { 5345a4b08SToby Isaac Vec diag; 6345a4b08SToby Isaac PetscBool diag_valid; 7345a4b08SToby Isaac Vec inv_diag; 8345a4b08SToby Isaac PetscBool inv_diag_valid; 9345a4b08SToby Isaac PetscObjectState diag_state, inv_diag_state; 10*c3e1b152SPierre Jolivet PetscInt *col; 11*c3e1b152SPierre Jolivet PetscScalar *val; 12345a4b08SToby Isaac } Mat_Diagonal; 13345a4b08SToby Isaac 14345a4b08SToby Isaac static PetscErrorCode MatDiagonalSetUpDiagonal(Mat A) 15345a4b08SToby Isaac { 16345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 17345a4b08SToby Isaac 18345a4b08SToby Isaac PetscFunctionBegin; 19345a4b08SToby Isaac if (!ctx->diag_valid) { 20345a4b08SToby Isaac PetscAssert(ctx->inv_diag_valid, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Neither diagonal nor inverse diagonal is in a valid state"); 21345a4b08SToby Isaac PetscCall(VecCopy(ctx->inv_diag, ctx->diag)); 22345a4b08SToby Isaac PetscCall(VecReciprocal(ctx->diag)); 23345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 24345a4b08SToby Isaac } 25345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 26345a4b08SToby Isaac } 27345a4b08SToby Isaac 28345a4b08SToby Isaac static PetscErrorCode MatDiagonalSetUpInverseDiagonal(Mat A) 29345a4b08SToby Isaac { 30345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 31345a4b08SToby Isaac 32345a4b08SToby Isaac PetscFunctionBegin; 33345a4b08SToby Isaac if (!ctx->inv_diag_valid) { 34345a4b08SToby Isaac PetscAssert(ctx->diag_valid, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Neither diagonal nor inverse diagonal is in a valid state"); 35345a4b08SToby Isaac PetscCall(VecCopy(ctx->diag, ctx->inv_diag)); 36345a4b08SToby Isaac PetscCall(VecReciprocal(ctx->inv_diag)); 37345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_TRUE; 38345a4b08SToby Isaac } 39345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 40345a4b08SToby Isaac } 41345a4b08SToby Isaac 42345a4b08SToby Isaac static PetscErrorCode MatAXPY_Diagonal(Mat Y, PetscScalar a, Mat X, MatStructure str) 43345a4b08SToby Isaac { 44345a4b08SToby Isaac Mat_Diagonal *yctx = (Mat_Diagonal *)Y->data; 45345a4b08SToby Isaac Mat_Diagonal *xctx = (Mat_Diagonal *)X->data; 46345a4b08SToby Isaac 47345a4b08SToby Isaac PetscFunctionBegin; 48345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(Y)); 49345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(X)); 50345a4b08SToby Isaac PetscCall(VecAXPY(yctx->diag, a, xctx->diag)); 51345a4b08SToby Isaac yctx->inv_diag_valid = PETSC_FALSE; 52345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 53345a4b08SToby Isaac } 54345a4b08SToby Isaac 55*c3e1b152SPierre Jolivet static PetscErrorCode MatGetRow_Diagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **vals) 56*c3e1b152SPierre Jolivet { 57*c3e1b152SPierre Jolivet Mat_Diagonal *mat = (Mat_Diagonal *)A->data; 58*c3e1b152SPierre Jolivet 59*c3e1b152SPierre Jolivet PetscFunctionBegin; 60*c3e1b152SPierre Jolivet if (ncols) *ncols = 1; 61*c3e1b152SPierre Jolivet if (cols) { 62*c3e1b152SPierre Jolivet if (!mat->col) PetscCall(PetscMalloc1(1, &mat->col)); 63*c3e1b152SPierre Jolivet *mat->col = row; 64*c3e1b152SPierre Jolivet *cols = mat->col; 65*c3e1b152SPierre Jolivet } 66*c3e1b152SPierre Jolivet if (vals) { 67*c3e1b152SPierre Jolivet const PetscScalar *v; 68*c3e1b152SPierre Jolivet 69*c3e1b152SPierre Jolivet if (!mat->val) PetscCall(PetscMalloc1(1, &mat->val)); 70*c3e1b152SPierre Jolivet PetscCall(VecGetArrayRead(mat->diag, &v)); 71*c3e1b152SPierre Jolivet *mat->val = v[row]; 72*c3e1b152SPierre Jolivet *vals = mat->val; 73*c3e1b152SPierre Jolivet PetscCall(VecRestoreArrayRead(mat->diag, &v)); 74*c3e1b152SPierre Jolivet } 75*c3e1b152SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 76*c3e1b152SPierre Jolivet } 77*c3e1b152SPierre Jolivet 78*c3e1b152SPierre Jolivet static PetscErrorCode MatRestoreRow_Diagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **vals) 79*c3e1b152SPierre Jolivet { 80*c3e1b152SPierre Jolivet PetscFunctionBegin; 81*c3e1b152SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 82*c3e1b152SPierre Jolivet } 83*c3e1b152SPierre Jolivet 84345a4b08SToby Isaac static PetscErrorCode MatMult_Diagonal(Mat A, Vec x, Vec y) 85345a4b08SToby Isaac { 86345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 87345a4b08SToby Isaac 88345a4b08SToby Isaac PetscFunctionBegin; 89345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(A)); 90345a4b08SToby Isaac PetscCall(VecPointwiseMult(y, ctx->diag, x)); 91345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 92345a4b08SToby Isaac } 93345a4b08SToby Isaac 94345a4b08SToby Isaac static PetscErrorCode MatMultAdd_Diagonal(Mat mat, Vec v1, Vec v2, Vec v3) 95345a4b08SToby Isaac { 96345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)mat->data; 97345a4b08SToby Isaac 98345a4b08SToby Isaac PetscFunctionBegin; 99345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(mat)); 100345a4b08SToby Isaac if (v2 != v3) { 101345a4b08SToby Isaac PetscCall(VecPointwiseMult(v3, ctx->diag, v1)); 102345a4b08SToby Isaac PetscCall(VecAXPY(v3, 1.0, v2)); 103345a4b08SToby Isaac } else { 104345a4b08SToby Isaac Vec w; 105345a4b08SToby Isaac PetscCall(VecDuplicate(v3, &w)); 106345a4b08SToby Isaac PetscCall(VecPointwiseMult(w, ctx->diag, v1)); 107345a4b08SToby Isaac PetscCall(VecAXPY(v3, 1.0, w)); 108345a4b08SToby Isaac PetscCall(VecDestroy(&w)); 109345a4b08SToby Isaac } 110345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 111345a4b08SToby Isaac } 112345a4b08SToby Isaac 113345a4b08SToby Isaac static PetscErrorCode MatNorm_Diagonal(Mat A, NormType type, PetscReal *nrm) 114345a4b08SToby Isaac { 115345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 116345a4b08SToby Isaac 117345a4b08SToby Isaac PetscFunctionBegin; 118345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(A)); 119345a4b08SToby Isaac type = (type == NORM_FROBENIUS) ? NORM_2 : NORM_INFINITY; 120345a4b08SToby Isaac PetscCall(VecNorm(ctx->diag, type, nrm)); 121345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 122345a4b08SToby Isaac } 123345a4b08SToby Isaac 124345a4b08SToby Isaac static PetscErrorCode MatDuplicate_Diagonal(Mat A, MatDuplicateOption op, Mat *B) 125345a4b08SToby Isaac { 126345a4b08SToby Isaac Mat_Diagonal *actx = (Mat_Diagonal *)A->data; 127345a4b08SToby Isaac 128345a4b08SToby Isaac PetscFunctionBegin; 129345a4b08SToby Isaac PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B)); 130345a4b08SToby Isaac PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 131345a4b08SToby Isaac PetscCall(MatSetBlockSizesFromMats(*B, A, A)); 132345a4b08SToby Isaac PetscCall(MatSetType(*B, MATDIAGONAL)); 133345a4b08SToby Isaac PetscCall(PetscLayoutReference(A->rmap, &(*B)->rmap)); 134345a4b08SToby Isaac PetscCall(PetscLayoutReference(A->cmap, &(*B)->cmap)); 135345a4b08SToby Isaac if (op == MAT_COPY_VALUES) { 136345a4b08SToby Isaac Mat_Diagonal *bctx = (Mat_Diagonal *)(*B)->data; 137345a4b08SToby Isaac 138345a4b08SToby Isaac PetscCall(MatSetUp(A)); 139345a4b08SToby Isaac PetscCall(MatSetUp(*B)); 140345a4b08SToby Isaac PetscCall(VecCopy(actx->diag, bctx->diag)); 141345a4b08SToby Isaac PetscCall(VecCopy(actx->inv_diag, bctx->inv_diag)); 142345a4b08SToby Isaac bctx->diag_valid = actx->diag_valid; 143345a4b08SToby Isaac bctx->inv_diag_valid = actx->inv_diag_valid; 144345a4b08SToby Isaac } 145345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 146345a4b08SToby Isaac } 147345a4b08SToby Isaac 148345a4b08SToby Isaac /*@ 149345a4b08SToby Isaac MatDiagonalGetDiagonal - Get the diagonal of a `MATDIAGONAL` 150345a4b08SToby Isaac 151345a4b08SToby Isaac Input Parameter: 152345a4b08SToby Isaac . A - the `MATDIAGONAL` 153345a4b08SToby Isaac 154345a4b08SToby Isaac Output Parameter: 155345a4b08SToby Isaac . diag - the `Vec` that defines the diagonal 156345a4b08SToby Isaac 157345a4b08SToby Isaac Level: developer 158345a4b08SToby Isaac 159345a4b08SToby Isaac Note: 160345a4b08SToby Isaac The user must call 161345a4b08SToby Isaac `MatDiagonalRestoreDiagonal()` before using the matrix again. 162345a4b08SToby Isaac 163345a4b08SToby Isaac For a copy of the diagonal values, rather than a reference, use `MatGetDiagonal()` 164345a4b08SToby Isaac 165345a4b08SToby 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()` 166345a4b08SToby Isaac 167be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()`, `MatGetDiagonal()` 168345a4b08SToby Isaac @*/ 169345a4b08SToby Isaac PetscErrorCode MatDiagonalGetDiagonal(Mat A, Vec *diag) 170345a4b08SToby Isaac { 171345a4b08SToby Isaac PetscFunctionBegin; 172345a4b08SToby Isaac PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 1734f572ea9SToby Isaac PetscAssertPointer(diag, 2); 174345a4b08SToby Isaac *diag = NULL; 175345a4b08SToby Isaac PetscUseMethod((PetscObject)A, "MatDiagonalGetDiagonal_C", (Mat, Vec *), (A, diag)); 176345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 177345a4b08SToby Isaac } 178345a4b08SToby Isaac 179345a4b08SToby Isaac static PetscErrorCode MatDiagonalGetDiagonal_Diagonal(Mat A, Vec *diag) 180345a4b08SToby Isaac { 181345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 182345a4b08SToby Isaac 183345a4b08SToby Isaac PetscFunctionBegin; 184345a4b08SToby Isaac PetscCall(MatSetUp(A)); 185345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(A)); 186345a4b08SToby Isaac *diag = ctx->diag; 187345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)*diag, &ctx->diag_state)); 188345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 189345a4b08SToby Isaac } 190345a4b08SToby Isaac 191345a4b08SToby Isaac /*@ 192345a4b08SToby Isaac MatDiagonalRestoreDiagonal - Restore the diagonal of a `MATDIAGONAL` 193345a4b08SToby Isaac 194345a4b08SToby Isaac Input Parameters: 195345a4b08SToby Isaac + A - the `MATDIAGONAL` 196345a4b08SToby Isaac - diag - the `Vec` obtained from `MatDiagonalGetDiagonal()` 197345a4b08SToby Isaac 198345a4b08SToby Isaac Level: developer 199345a4b08SToby Isaac 200345a4b08SToby Isaac Note: 201345a4b08SToby Isaac Use `MatDiagonalSet()` to change the values by copy, rather than reference. 202345a4b08SToby Isaac 203be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalGetDiagonal()` 204345a4b08SToby Isaac @*/ 205345a4b08SToby Isaac PetscErrorCode MatDiagonalRestoreDiagonal(Mat A, Vec *diag) 206345a4b08SToby Isaac { 207345a4b08SToby Isaac PetscFunctionBegin; 208345a4b08SToby Isaac PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2094f572ea9SToby Isaac PetscAssertPointer(diag, 2); 210345a4b08SToby Isaac PetscUseMethod((PetscObject)A, "MatDiagonalRestoreDiagonal_C", (Mat, Vec *), (A, diag)); 211345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 212345a4b08SToby Isaac } 213345a4b08SToby Isaac 214345a4b08SToby Isaac static PetscErrorCode MatDiagonalRestoreDiagonal_Diagonal(Mat A, Vec *diag) 215345a4b08SToby Isaac { 216345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 217345a4b08SToby Isaac PetscObjectState diag_state; 218345a4b08SToby Isaac 219345a4b08SToby Isaac PetscFunctionBegin; 220345a4b08SToby Isaac PetscCheck(ctx->diag == *diag, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Restored a different diagonal vector"); 221345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 222345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)*diag, &diag_state)); 223345a4b08SToby Isaac if (ctx->diag_state != diag_state) { 224345a4b08SToby Isaac PetscCall(PetscObjectStateIncrease((PetscObject)A)); 225345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 226345a4b08SToby Isaac } 227345a4b08SToby Isaac *diag = NULL; 228345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 229345a4b08SToby Isaac } 230345a4b08SToby Isaac 231345a4b08SToby Isaac /*@ 232345a4b08SToby Isaac MatDiagonalGetInverseDiagonal - Get the inverse diagonal of a `MATDIAGONAL` 233345a4b08SToby Isaac 234345a4b08SToby Isaac Input Parameter: 235345a4b08SToby Isaac . A - the `MATDIAGONAL` 236345a4b08SToby Isaac 237345a4b08SToby Isaac Output Parameter: 238345a4b08SToby Isaac . inv_diag - the `Vec` that defines the inverse diagonal 239345a4b08SToby Isaac 240345a4b08SToby Isaac Level: developer 241345a4b08SToby Isaac 242345a4b08SToby Isaac Note: 243345a4b08SToby Isaac The user must call 244345a4b08SToby Isaac `MatDiagonalRestoreInverseDiagonal()` before using the matrix again. 245345a4b08SToby Isaac 246345a4b08SToby Isaac If a matrix is created only to call `MatSolve()` (which happens for `MATLMVMDIAGBROYDEN`), 247345a4b08SToby Isaac using `MatDiagonalGetInverseDiagonal()` and `MatDiagonalRestoreInverseDiagonal()` avoids copies 248345a4b08SToby Isaac and avoids any call to `VecReciprocal()`. 249345a4b08SToby Isaac 250be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MATLMVMBROYDEN`, `MatSolve()` 251345a4b08SToby Isaac @*/ 252345a4b08SToby Isaac PetscErrorCode MatDiagonalGetInverseDiagonal(Mat A, Vec *inv_diag) 253345a4b08SToby Isaac { 254345a4b08SToby Isaac PetscFunctionBegin; 255345a4b08SToby Isaac PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2564f572ea9SToby Isaac PetscAssertPointer(inv_diag, 2); 257345a4b08SToby Isaac *inv_diag = NULL; 258345a4b08SToby Isaac PetscUseMethod((PetscObject)A, "MatDiagonalGetInverseDiagonal_C", (Mat, Vec *), (A, inv_diag)); 259345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 260345a4b08SToby Isaac } 261345a4b08SToby Isaac 262345a4b08SToby Isaac static PetscErrorCode MatDiagonalGetInverseDiagonal_Diagonal(Mat A, Vec *inv_diag) 263345a4b08SToby Isaac { 264345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 265345a4b08SToby Isaac 266345a4b08SToby Isaac PetscFunctionBegin; 267345a4b08SToby Isaac PetscCall(MatSetUp(A)); 268345a4b08SToby Isaac PetscCall(MatDiagonalSetUpInverseDiagonal(A)); 269345a4b08SToby Isaac *inv_diag = ctx->inv_diag; 270345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)*inv_diag, &ctx->inv_diag_state)); 271345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 272345a4b08SToby Isaac } 273345a4b08SToby Isaac 274345a4b08SToby Isaac /*@ 275345a4b08SToby Isaac MatDiagonalRestoreInverseDiagonal - Restore the inverse diagonal of a `MATDIAGONAL` 276345a4b08SToby Isaac 277345a4b08SToby Isaac Input Parameters: 278345a4b08SToby Isaac + A - the `MATDIAGONAL` 279345a4b08SToby Isaac - inv_diag - the `Vec` obtained from `MatDiagonalGetInverseDiagonal()` 280345a4b08SToby Isaac 281345a4b08SToby Isaac Level: developer 282345a4b08SToby Isaac 283be50c303SSatish Balay .seealso: [](ch_matrices), `MATDIAGONAL`, `MatCreateDiagonal()`, `MatDiagonalGetInverseDiagonal()` 284345a4b08SToby Isaac @*/ 285345a4b08SToby Isaac PetscErrorCode MatDiagonalRestoreInverseDiagonal(Mat A, Vec *inv_diag) 286345a4b08SToby Isaac { 287345a4b08SToby Isaac PetscFunctionBegin; 288345a4b08SToby Isaac PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2894f572ea9SToby Isaac PetscAssertPointer(inv_diag, 2); 290345a4b08SToby Isaac PetscUseMethod((PetscObject)A, "MatDiagonalRestoreInverseDiagonal_C", (Mat, Vec *), (A, inv_diag)); 291345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 292345a4b08SToby Isaac } 293345a4b08SToby Isaac 294345a4b08SToby Isaac static PetscErrorCode MatDiagonalRestoreInverseDiagonal_Diagonal(Mat A, Vec *inv_diag) 295345a4b08SToby Isaac { 296345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 297345a4b08SToby Isaac PetscObjectState inv_diag_state; 298345a4b08SToby Isaac 299345a4b08SToby Isaac PetscFunctionBegin; 300345a4b08SToby Isaac PetscCheck(ctx->inv_diag == *inv_diag, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Restored a different diagonal vector"); 301345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_TRUE; 302345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)*inv_diag, &inv_diag_state)); 303345a4b08SToby Isaac if (ctx->inv_diag_state != inv_diag_state) { 304345a4b08SToby Isaac PetscCall(PetscObjectStateIncrease((PetscObject)A)); 305345a4b08SToby Isaac ctx->diag_valid = PETSC_FALSE; 306345a4b08SToby Isaac } 307345a4b08SToby Isaac *inv_diag = NULL; 308345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 309345a4b08SToby Isaac } 310345a4b08SToby Isaac 311345a4b08SToby Isaac static PetscErrorCode MatDestroy_Diagonal(Mat mat) 312345a4b08SToby Isaac { 313345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)mat->data; 314345a4b08SToby Isaac 315345a4b08SToby Isaac PetscFunctionBegin; 316345a4b08SToby Isaac PetscCall(VecDestroy(&ctx->diag)); 317345a4b08SToby Isaac PetscCall(VecDestroy(&ctx->inv_diag)); 318*c3e1b152SPierre Jolivet PetscCall(PetscFree(ctx->col)); 319*c3e1b152SPierre Jolivet PetscCall(PetscFree(ctx->val)); 320345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalGetDiagonal_C", NULL)); 321345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalRestoreDiagonal_C", NULL)); 322345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalGetInverseDiagonal_C", NULL)); 323345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalRestoreInverseDiagonal_C", NULL)); 324345a4b08SToby Isaac PetscCall(PetscFree(mat->data)); 325345a4b08SToby Isaac mat->structural_symmetry_eternal = PETSC_FALSE; 326345a4b08SToby Isaac mat->symmetry_eternal = PETSC_FALSE; 327345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 328345a4b08SToby Isaac } 329345a4b08SToby Isaac 330345a4b08SToby Isaac static PetscErrorCode MatView_Diagonal(Mat J, PetscViewer viewer) 331345a4b08SToby Isaac { 332345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)J->data; 333345a4b08SToby Isaac PetscBool iascii; 334345a4b08SToby Isaac 335345a4b08SToby Isaac PetscFunctionBegin; 336345a4b08SToby Isaac PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 337345a4b08SToby Isaac if (iascii) { 338345a4b08SToby Isaac PetscViewerFormat format; 339345a4b08SToby Isaac 340345a4b08SToby Isaac PetscCall(PetscViewerGetFormat(viewer, &format)); 341345a4b08SToby Isaac if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS); 342345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 343345a4b08SToby Isaac PetscCall(VecView(ctx->diag, viewer)); 344345a4b08SToby Isaac } 345345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 346345a4b08SToby Isaac } 347345a4b08SToby Isaac 348345a4b08SToby Isaac static PetscErrorCode MatGetDiagonal_Diagonal(Mat J, Vec x) 349345a4b08SToby Isaac { 350345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)J->data; 351345a4b08SToby Isaac 352345a4b08SToby Isaac PetscFunctionBegin; 353345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 354345a4b08SToby Isaac PetscCall(VecCopy(ctx->diag, x)); 355345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 356345a4b08SToby Isaac } 357345a4b08SToby Isaac 358345a4b08SToby Isaac static PetscErrorCode MatDiagonalSet_Diagonal(Mat J, Vec D, InsertMode is) 359345a4b08SToby Isaac { 360345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)J->data; 361345a4b08SToby Isaac 362345a4b08SToby Isaac PetscFunctionBegin; 363345a4b08SToby Isaac switch (is) { 364345a4b08SToby Isaac case ADD_VALUES: 365345a4b08SToby Isaac case ADD_ALL_VALUES: 366345a4b08SToby Isaac case ADD_BC_VALUES: 367345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 368345a4b08SToby Isaac PetscCall(VecAXPY(ctx->diag, 1.0, D)); 369345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 370345a4b08SToby Isaac break; 371345a4b08SToby Isaac case INSERT_VALUES: 372345a4b08SToby Isaac case INSERT_BC_VALUES: 373345a4b08SToby Isaac case INSERT_ALL_VALUES: 374345a4b08SToby Isaac PetscCall(MatSetUp(J)); 375345a4b08SToby Isaac PetscCall(VecCopy(D, ctx->diag)); 376345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 377345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 378345a4b08SToby Isaac break; 379345a4b08SToby Isaac case MAX_VALUES: 380345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 381345a4b08SToby Isaac PetscCall(VecPointwiseMax(ctx->diag, D, ctx->diag)); 382345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 383345a4b08SToby Isaac break; 384345a4b08SToby Isaac case MIN_VALUES: 385345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(J)); 386345a4b08SToby Isaac PetscCall(VecPointwiseMin(ctx->diag, D, ctx->diag)); 387345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 388345a4b08SToby Isaac break; 389345a4b08SToby Isaac case NOT_SET_VALUES: 390345a4b08SToby Isaac break; 391345a4b08SToby Isaac } 392345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 393345a4b08SToby Isaac } 394345a4b08SToby Isaac 395345a4b08SToby Isaac static PetscErrorCode MatShift_Diagonal(Mat Y, PetscScalar a) 396345a4b08SToby Isaac { 397345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data; 398345a4b08SToby Isaac 399345a4b08SToby Isaac PetscFunctionBegin; 400345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(Y)); 401345a4b08SToby Isaac PetscCall(VecShift(ctx->diag, a)); 402345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 403345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 404345a4b08SToby Isaac } 405345a4b08SToby Isaac 406345a4b08SToby Isaac static PetscErrorCode MatScale_Diagonal(Mat Y, PetscScalar a) 407345a4b08SToby Isaac { 408345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data; 409345a4b08SToby Isaac 410345a4b08SToby Isaac PetscFunctionBegin; 411345a4b08SToby Isaac PetscCall(MatSetUp(Y)); 412345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(Y)); 413345a4b08SToby Isaac PetscCall(VecScale(ctx->diag, a)); 414345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 415345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 416345a4b08SToby Isaac } 417345a4b08SToby Isaac 418345a4b08SToby Isaac static PetscErrorCode MatDiagonalScale_Diagonal(Mat Y, Vec l, Vec r) 419345a4b08SToby Isaac { 420345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data; 421345a4b08SToby Isaac 422345a4b08SToby Isaac PetscFunctionBegin; 423345a4b08SToby Isaac PetscCall(MatDiagonalSetUpDiagonal(Y)); 424345a4b08SToby Isaac if (l) { 425345a4b08SToby Isaac PetscCall(VecPointwiseMult(ctx->diag, ctx->diag, l)); 426345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 427345a4b08SToby Isaac } 428345a4b08SToby Isaac if (r) { 429345a4b08SToby Isaac PetscCall(VecPointwiseMult(ctx->diag, ctx->diag, r)); 430345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 431345a4b08SToby Isaac } 432345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 433345a4b08SToby Isaac } 434345a4b08SToby Isaac 435345a4b08SToby Isaac static PetscErrorCode MatSetUp_Diagonal(Mat A) 436345a4b08SToby Isaac { 437345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)A->data; 438345a4b08SToby Isaac 439345a4b08SToby Isaac PetscFunctionBegin; 440345a4b08SToby Isaac if (!ctx->diag) { 441345a4b08SToby Isaac PetscCall(PetscLayoutSetUp(A->cmap)); 442345a4b08SToby Isaac PetscCall(PetscLayoutSetUp(A->rmap)); 443345a4b08SToby Isaac PetscCall(MatCreateVecs(A, &ctx->diag, NULL)); 444345a4b08SToby Isaac PetscCall(VecDuplicate(ctx->diag, &ctx->inv_diag)); 445345a4b08SToby Isaac PetscCall(VecZeroEntries(ctx->diag)); 446345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 447345a4b08SToby Isaac } 448345a4b08SToby Isaac A->assembled = PETSC_TRUE; 449345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 450345a4b08SToby Isaac } 451345a4b08SToby Isaac 452345a4b08SToby Isaac static PetscErrorCode MatZeroEntries_Diagonal(Mat Y) 453345a4b08SToby Isaac { 454345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)Y->data; 455345a4b08SToby Isaac 456345a4b08SToby Isaac PetscFunctionBegin; 457345a4b08SToby Isaac PetscCall(MatSetUp(Y)); 458345a4b08SToby Isaac PetscCall(VecZeroEntries(ctx->diag)); 459345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 460345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 461345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 462345a4b08SToby Isaac } 463345a4b08SToby Isaac 46466976f2fSJacob Faibussowitsch static PetscErrorCode MatSolve_Diagonal(Mat matin, Vec b, Vec x) 465345a4b08SToby Isaac { 466345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)matin->data; 467345a4b08SToby Isaac 468345a4b08SToby Isaac PetscFunctionBegin; 469345a4b08SToby Isaac PetscCall(MatDiagonalSetUpInverseDiagonal(matin)); 470345a4b08SToby Isaac PetscCall(VecPointwiseMult(x, b, ctx->inv_diag)); 471345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 472345a4b08SToby Isaac } 473345a4b08SToby Isaac 47466976f2fSJacob Faibussowitsch static PetscErrorCode MatGetInfo_Diagonal(Mat A, MatInfoType flag, MatInfo *info) 475345a4b08SToby Isaac { 476345a4b08SToby Isaac PetscFunctionBegin; 477345a4b08SToby Isaac info->block_size = 1.0; 478345a4b08SToby Isaac info->nz_allocated = A->cmap->N; 479345a4b08SToby Isaac info->nz_used = A->cmap->N; 480345a4b08SToby Isaac info->nz_unneeded = 0.0; 481345a4b08SToby Isaac info->assemblies = A->num_ass; 482345a4b08SToby Isaac info->mallocs = 0.0; 483345a4b08SToby Isaac info->memory = 0; 484345a4b08SToby Isaac info->fill_ratio_given = 0; 485345a4b08SToby Isaac info->fill_ratio_needed = 0; 486345a4b08SToby Isaac info->factor_mallocs = 0; 487345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 488345a4b08SToby Isaac } 489345a4b08SToby Isaac 490345a4b08SToby Isaac /*@ 491345a4b08SToby Isaac MatCreateDiagonal - Creates a matrix defined by a given vector along its diagonal. 492345a4b08SToby Isaac 493345a4b08SToby Isaac Collective 494345a4b08SToby Isaac 495345a4b08SToby Isaac Input Parameter: 496345a4b08SToby Isaac . diag - vector for the diagonal 497345a4b08SToby Isaac 498345a4b08SToby Isaac Output Parameter: 499345a4b08SToby Isaac . J - the diagonal matrix 500345a4b08SToby Isaac 501345a4b08SToby Isaac Level: advanced 502345a4b08SToby Isaac 503345a4b08SToby Isaac Notes: 504345a4b08SToby Isaac Only supports square matrices with the same number of local rows and columns. 505345a4b08SToby Isaac 506345a4b08SToby Isaac The input vector `diag` will be referenced internally: any changes to `diag` 507345a4b08SToby Isaac will affect the matrix `J`. 508345a4b08SToby Isaac 509be50c303SSatish Balay .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatSolve()` 510345a4b08SToby Isaac `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()` 511345a4b08SToby Isaac @*/ 512345a4b08SToby Isaac PetscErrorCode MatCreateDiagonal(Vec diag, Mat *J) 513345a4b08SToby Isaac { 514345a4b08SToby Isaac PetscFunctionBegin; 515345a4b08SToby Isaac PetscValidHeaderSpecific(diag, VEC_CLASSID, 1); 516345a4b08SToby Isaac PetscCall(MatCreate(PetscObjectComm((PetscObject)diag), J)); 517345a4b08SToby Isaac PetscInt m, M; 518345a4b08SToby Isaac PetscCall(VecGetLocalSize(diag, &m)); 519345a4b08SToby Isaac PetscCall(VecGetSize(diag, &M)); 520345a4b08SToby Isaac PetscCall(MatSetSizes(*J, m, m, M, M)); 521345a4b08SToby Isaac PetscCall(MatSetType(*J, MATDIAGONAL)); 522345a4b08SToby Isaac 523345a4b08SToby Isaac PetscLayout map; 524345a4b08SToby Isaac PetscCall(VecGetLayout(diag, &map)); 525345a4b08SToby Isaac PetscCall(MatSetLayouts(*J, map, map)); 526345a4b08SToby Isaac Mat_Diagonal *ctx = (Mat_Diagonal *)(*J)->data; 527345a4b08SToby Isaac PetscCall(PetscObjectReference((PetscObject)diag)); 528345a4b08SToby Isaac PetscCall(VecDestroy(&ctx->diag)); 529345a4b08SToby Isaac PetscCall(VecDestroy(&ctx->inv_diag)); 530345a4b08SToby Isaac ctx->diag = diag; 531345a4b08SToby Isaac ctx->diag_valid = PETSC_TRUE; 532345a4b08SToby Isaac ctx->inv_diag_valid = PETSC_FALSE; 533345a4b08SToby Isaac VecType type; 534345a4b08SToby Isaac PetscCall(VecDuplicate(ctx->diag, &ctx->inv_diag)); 535345a4b08SToby Isaac PetscCall(VecGetType(diag, &type)); 536345a4b08SToby Isaac PetscCall(PetscFree((*J)->defaultvectype)); 537345a4b08SToby Isaac PetscCall(PetscStrallocpy(type, &(*J)->defaultvectype)); 538345a4b08SToby Isaac PetscCall(MatSetUp(*J)); 539*c3e1b152SPierre Jolivet ctx->col = NULL; 540*c3e1b152SPierre Jolivet ctx->val = NULL; 541345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 542345a4b08SToby Isaac } 543345a4b08SToby Isaac 544345a4b08SToby Isaac /*MC 545345a4b08SToby Isaac MATDIAGONAL - MATDIAGONAL = "diagonal" - A diagonal matrix type with the diagonal implemented as a `Vec`. Useful for 546345a4b08SToby Isaac cases where `VecPointwiseMult()` or `VecPointwiseDivide()` should be thought of as the actions of a linear operator. 547345a4b08SToby Isaac 548345a4b08SToby Isaac Level: advanced 549345a4b08SToby Isaac 550be50c303SSatish Balay .seealso: [](ch_matrices), `Mat`, `MatCreateDiagonal()`, `MatDiagonalRestoreInverseDiagonal()`, `MatDiagonalGetDiagonal()`, `MatDiagonalRestoreDiagonal()`, `MatDiagonalGetInverseDiagonal()` 551345a4b08SToby Isaac M*/ 552345a4b08SToby Isaac PETSC_INTERN PetscErrorCode MatCreate_Diagonal(Mat A) 553345a4b08SToby Isaac { 554345a4b08SToby Isaac Mat_Diagonal *ctx; 555345a4b08SToby Isaac 556345a4b08SToby Isaac PetscFunctionBegin; 557345a4b08SToby Isaac PetscCall(PetscNew(&ctx)); 558345a4b08SToby Isaac A->data = (void *)ctx; 559345a4b08SToby Isaac 560345a4b08SToby Isaac A->structurally_symmetric = PETSC_BOOL3_TRUE; 561345a4b08SToby Isaac A->structural_symmetry_eternal = PETSC_TRUE; 562345a4b08SToby Isaac A->symmetry_eternal = PETSC_TRUE; 563345a4b08SToby Isaac A->symmetric = PETSC_BOOL3_TRUE; 564345a4b08SToby Isaac if (!PetscDefined(USE_COMPLEX)) A->hermitian = PETSC_BOOL3_TRUE; 565345a4b08SToby Isaac 566*c3e1b152SPierre Jolivet A->ops->getrow = MatGetRow_Diagonal; 567*c3e1b152SPierre Jolivet A->ops->restorerow = MatRestoreRow_Diagonal; 568345a4b08SToby Isaac A->ops->mult = MatMult_Diagonal; 569345a4b08SToby Isaac A->ops->multadd = MatMultAdd_Diagonal; 570345a4b08SToby Isaac A->ops->multtranspose = MatMult_Diagonal; 571345a4b08SToby Isaac A->ops->multtransposeadd = MatMultAdd_Diagonal; 572345a4b08SToby Isaac A->ops->norm = MatNorm_Diagonal; 573345a4b08SToby Isaac A->ops->duplicate = MatDuplicate_Diagonal; 574345a4b08SToby Isaac A->ops->solve = MatSolve_Diagonal; 575345a4b08SToby Isaac A->ops->solvetranspose = MatSolve_Diagonal; 576345a4b08SToby Isaac A->ops->shift = MatShift_Diagonal; 577345a4b08SToby Isaac A->ops->scale = MatScale_Diagonal; 578345a4b08SToby Isaac A->ops->diagonalscale = MatDiagonalScale_Diagonal; 579345a4b08SToby Isaac A->ops->getdiagonal = MatGetDiagonal_Diagonal; 580345a4b08SToby Isaac A->ops->diagonalset = MatDiagonalSet_Diagonal; 581345a4b08SToby Isaac A->ops->view = MatView_Diagonal; 582345a4b08SToby Isaac A->ops->zeroentries = MatZeroEntries_Diagonal; 583345a4b08SToby Isaac A->ops->destroy = MatDestroy_Diagonal; 584345a4b08SToby Isaac A->ops->getinfo = MatGetInfo_Diagonal; 585345a4b08SToby Isaac A->ops->axpy = MatAXPY_Diagonal; 586345a4b08SToby Isaac A->ops->setup = MatSetUp_Diagonal; 587345a4b08SToby Isaac 588345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalGetDiagonal_C", MatDiagonalGetDiagonal_Diagonal)); 589345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalRestoreDiagonal_C", MatDiagonalRestoreDiagonal_Diagonal)); 590345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalGetInverseDiagonal_C", MatDiagonalGetInverseDiagonal_Diagonal)); 591345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatDiagonalRestoreInverseDiagonal_C", MatDiagonalRestoreInverseDiagonal_Diagonal)); 592345a4b08SToby Isaac PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATDIAGONAL)); 593345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 594345a4b08SToby Isaac } 595