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