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