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