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