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