1 /* 2 Defines a preconditioner that can consist of a collection of PCs 3 */ 4 #include <petsc/private/pcimpl.h> 5 #include <petscksp.h> /*I "petscksp.h" I*/ 6 7 typedef struct _PC_CompositeLink *PC_CompositeLink; 8 struct _PC_CompositeLink { 9 PC pc; 10 PC_CompositeLink next; 11 PC_CompositeLink previous; 12 }; 13 14 typedef struct { 15 PC_CompositeLink head; 16 PCCompositeType type; 17 Vec work1; 18 Vec work2; 19 PetscScalar alpha; 20 Mat alpha_mat; 21 } PC_Composite; 22 23 static PetscErrorCode PCApply_Composite_Multiplicative(PC pc, Vec x, Vec y) 24 { 25 PC_Composite *jac = (PC_Composite *)pc->data; 26 PC_CompositeLink next = jac->head; 27 Mat mat = pc->pmat; 28 29 PetscFunctionBegin; 30 PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 31 32 /* Set the reuse flag on children PCs */ 33 while (next) { 34 PetscCall(PCSetReusePreconditioner(next->pc, pc->reusepreconditioner)); 35 next = next->next; 36 } 37 next = jac->head; 38 39 if (next->next && !jac->work2) { /* allocate second work vector */ 40 PetscCall(VecDuplicate(jac->work1, &jac->work2)); 41 } 42 if (pc->useAmat) mat = pc->mat; 43 PetscCall(PCApply(next->pc, x, y)); /* y <- B x */ 44 while (next->next) { 45 next = next->next; 46 PetscCall(MatMult(mat, y, jac->work1)); /* work1 <- A y */ 47 PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x)); /* work2 <- x - work1 */ 48 PetscCall(PCApply(next->pc, jac->work2, jac->work1)); /* work1 <- C work2 */ 49 PetscCall(VecAXPY(y, 1.0, jac->work1)); /* y <- y + work1 = B x + C (x - A B x) = (B + C (1 - A B)) x */ 50 } 51 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 52 while (next->previous) { 53 next = next->previous; 54 PetscCall(MatMult(mat, y, jac->work1)); 55 PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x)); 56 PetscCall(PCApply(next->pc, jac->work2, jac->work1)); 57 PetscCall(VecAXPY(y, 1.0, jac->work1)); 58 } 59 } 60 PetscFunctionReturn(PETSC_SUCCESS); 61 } 62 63 static PetscErrorCode PCApplyTranspose_Composite_Multiplicative(PC pc, Vec x, Vec y) 64 { 65 PC_Composite *jac = (PC_Composite *)pc->data; 66 PC_CompositeLink next = jac->head; 67 Mat mat = pc->pmat; 68 69 PetscFunctionBegin; 70 PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 71 if (next->next && !jac->work2) { /* allocate second work vector */ 72 PetscCall(VecDuplicate(jac->work1, &jac->work2)); 73 } 74 if (pc->useAmat) mat = pc->mat; 75 /* locate last PC */ 76 while (next->next) next = next->next; 77 PetscCall(PCApplyTranspose(next->pc, x, y)); 78 while (next->previous) { 79 next = next->previous; 80 PetscCall(MatMultTranspose(mat, y, jac->work1)); 81 PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x)); 82 PetscCall(PCApplyTranspose(next->pc, jac->work2, jac->work1)); 83 PetscCall(VecAXPY(y, 1.0, jac->work1)); 84 } 85 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 86 next = jac->head; 87 while (next->next) { 88 next = next->next; 89 PetscCall(MatMultTranspose(mat, y, jac->work1)); 90 PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x)); 91 PetscCall(PCApplyTranspose(next->pc, jac->work2, jac->work1)); 92 PetscCall(VecAXPY(y, 1.0, jac->work1)); 93 } 94 } 95 PetscFunctionReturn(PETSC_SUCCESS); 96 } 97 98 /* 99 This is very special for a matrix of the form alpha I + R + S 100 where first preconditioner is built from alpha I + S and second from 101 alpha I + R 102 */ 103 static PetscErrorCode PCApply_Composite_Special(PC pc, Vec x, Vec y) 104 { 105 PC_Composite *jac = (PC_Composite *)pc->data; 106 PC_CompositeLink next = jac->head; 107 108 PetscFunctionBegin; 109 PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 110 PetscCheck(next->next && !next->next->next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Special composite preconditioners requires exactly two PCs"); 111 112 /* Set the reuse flag on children PCs */ 113 PetscCall(PCSetReusePreconditioner(next->pc, pc->reusepreconditioner)); 114 PetscCall(PCSetReusePreconditioner(next->next->pc, pc->reusepreconditioner)); 115 116 PetscCall(PCApply(next->pc, x, jac->work1)); 117 if (jac->alpha_mat) { 118 if (!jac->work2) PetscCall(VecDuplicate(jac->work1, &jac->work2)); 119 PetscCall(MatMult(jac->alpha_mat, jac->work1, jac->work2)); 120 PetscCall(PCApply(next->next->pc, jac->work2, y)); 121 } else PetscCall(PCApply(next->next->pc, jac->work1, y)); 122 PetscFunctionReturn(PETSC_SUCCESS); 123 } 124 125 static PetscErrorCode PCApply_Composite_Additive(PC pc, Vec x, Vec y) 126 { 127 PC_Composite *jac = (PC_Composite *)pc->data; 128 PC_CompositeLink next = jac->head; 129 130 PetscFunctionBegin; 131 PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 132 133 /* Set the reuse flag on children PCs */ 134 while (next) { 135 PetscCall(PCSetReusePreconditioner(next->pc, pc->reusepreconditioner)); 136 next = next->next; 137 } 138 next = jac->head; 139 140 PetscCall(PCApply(next->pc, x, y)); 141 while (next->next) { 142 next = next->next; 143 PetscCall(PCApply(next->pc, x, jac->work1)); 144 PetscCall(VecAXPY(y, 1.0, jac->work1)); 145 } 146 PetscFunctionReturn(PETSC_SUCCESS); 147 } 148 149 static PetscErrorCode PCApplyTranspose_Composite_Additive(PC pc, Vec x, Vec y) 150 { 151 PC_Composite *jac = (PC_Composite *)pc->data; 152 PC_CompositeLink next = jac->head; 153 154 PetscFunctionBegin; 155 PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 156 PetscCall(PCApplyTranspose(next->pc, x, y)); 157 while (next->next) { 158 next = next->next; 159 PetscCall(PCApplyTranspose(next->pc, x, jac->work1)); 160 PetscCall(VecAXPY(y, 1.0, jac->work1)); 161 } 162 PetscFunctionReturn(PETSC_SUCCESS); 163 } 164 165 static PetscErrorCode PCSetUp_Composite(PC pc) 166 { 167 PC_Composite *jac = (PC_Composite *)pc->data; 168 PC_CompositeLink next = jac->head; 169 DM dm; 170 171 PetscFunctionBegin; 172 if (!jac->work1) PetscCall(MatCreateVecs(pc->pmat, &jac->work1, NULL)); 173 PetscCall(PCGetDM(pc, &dm)); 174 while (next) { 175 if (!next->pc->dm) PetscCall(PCSetDM(next->pc, dm)); 176 if (!next->pc->mat) PetscCall(PCSetOperators(next->pc, pc->mat, pc->pmat)); 177 next = next->next; 178 } 179 PetscFunctionReturn(PETSC_SUCCESS); 180 } 181 182 static PetscErrorCode PCSetUpOnBlocks_Composite(PC pc) 183 { 184 PC_Composite *jac = (PC_Composite *)pc->data; 185 PC_CompositeLink next = jac->head; 186 PCFailedReason reason; 187 188 PetscFunctionBegin; 189 while (next) { 190 PetscCall(PCSetUp(next->pc)); 191 PetscCall(PCGetFailedReasonRank(next->pc, &reason)); 192 if (reason) pc->failedreason = reason; 193 next = next->next; 194 } 195 PetscFunctionReturn(PETSC_SUCCESS); 196 } 197 198 static PetscErrorCode PCReset_Composite(PC pc) 199 { 200 PC_Composite *jac = (PC_Composite *)pc->data; 201 PC_CompositeLink next = jac->head; 202 203 PetscFunctionBegin; 204 while (next) { 205 PetscCall(PCReset(next->pc)); 206 next = next->next; 207 } 208 PetscCall(VecDestroy(&jac->work1)); 209 PetscCall(VecDestroy(&jac->work2)); 210 PetscCall(MatDestroy(&jac->alpha_mat)); 211 PetscFunctionReturn(PETSC_SUCCESS); 212 } 213 214 static PetscErrorCode PCDestroy_Composite(PC pc) 215 { 216 PC_Composite *jac = (PC_Composite *)pc->data; 217 PC_CompositeLink next = jac->head, next_tmp; 218 219 PetscFunctionBegin; 220 PetscCall(PCReset_Composite(pc)); 221 while (next) { 222 PetscCall(PCDestroy(&next->pc)); 223 next_tmp = next; 224 next = next->next; 225 PetscCall(PetscFree(next_tmp)); 226 } 227 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSetType_C", NULL)); 228 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetType_C", NULL)); 229 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPCType_C", NULL)); 230 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPC_C", NULL)); 231 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetNumberPC_C", NULL)); 232 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetPC_C", NULL)); 233 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlpha_C", NULL)); 234 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlphaMat_C", NULL)); 235 PetscCall(PetscFree(pc->data)); 236 PetscFunctionReturn(PETSC_SUCCESS); 237 } 238 239 static PetscErrorCode PCSetFromOptions_Composite(PC pc, PetscOptionItems *PetscOptionsObject) 240 { 241 PC_Composite *jac = (PC_Composite *)pc->data; 242 PetscInt nmax = 8, i; 243 PC_CompositeLink next; 244 char *pcs[8]; 245 PetscBool flg; 246 247 PetscFunctionBegin; 248 PetscOptionsHeadBegin(PetscOptionsObject, "Composite preconditioner options"); 249 PetscCall(PetscOptionsEnum("-pc_composite_type", "Type of composition", "PCCompositeSetType", PCCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&jac->type, &flg)); 250 if (flg) PetscCall(PCCompositeSetType(pc, jac->type)); 251 PetscCall(PetscOptionsStringArray("-pc_composite_pcs", "List of composite solvers", "PCCompositeAddPCType", pcs, &nmax, &flg)); 252 if (flg) { 253 for (i = 0; i < nmax; i++) { 254 PetscCall(PCCompositeAddPCType(pc, pcs[i])); 255 PetscCall(PetscFree(pcs[i])); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */ 256 } 257 } 258 PetscOptionsHeadEnd(); 259 260 next = jac->head; 261 while (next) { 262 PetscCall(PCSetFromOptions(next->pc)); 263 next = next->next; 264 } 265 PetscFunctionReturn(PETSC_SUCCESS); 266 } 267 268 static PetscErrorCode PCView_Composite(PC pc, PetscViewer viewer) 269 { 270 PC_Composite *jac = (PC_Composite *)pc->data; 271 PC_CompositeLink next = jac->head; 272 PetscBool iascii; 273 274 PetscFunctionBegin; 275 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 276 if (iascii) { 277 PetscCall(PetscViewerASCIIPrintf(viewer, "Composite PC type - %s\n", PCCompositeTypes[jac->type])); 278 PetscCall(PetscViewerASCIIPrintf(viewer, "PCs on composite preconditioner follow\n")); 279 PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------\n")); 280 } 281 if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer)); 282 while (next) { 283 PetscCall(PCView(next->pc, viewer)); 284 next = next->next; 285 } 286 if (iascii) { 287 PetscCall(PetscViewerASCIIPopTab(viewer)); 288 PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------\n")); 289 } 290 PetscFunctionReturn(PETSC_SUCCESS); 291 } 292 293 static PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc, PetscScalar alpha) 294 { 295 PC_Composite *jac = (PC_Composite *)pc->data; 296 297 PetscFunctionBegin; 298 jac->alpha = alpha; 299 PetscFunctionReturn(PETSC_SUCCESS); 300 } 301 302 static PetscErrorCode PCCompositeSpecialSetAlphaMat_Composite(PC pc, Mat alpha_mat) 303 { 304 PC_Composite *jac = (PC_Composite *)pc->data; 305 306 PetscFunctionBegin; 307 if (alpha_mat) { 308 PetscValidHeaderSpecific(alpha_mat, MAT_CLASSID, 2); 309 PetscCall(PetscObjectReference((PetscObject)alpha_mat)); 310 } 311 PetscCall(MatDestroy(&jac->alpha_mat)); 312 jac->alpha_mat = alpha_mat; 313 PetscFunctionReturn(PETSC_SUCCESS); 314 } 315 316 static PetscErrorCode PCCompositeSetType_Composite(PC pc, PCCompositeType type) 317 { 318 PC_Composite *jac = (PC_Composite *)pc->data; 319 320 PetscFunctionBegin; 321 if (type == PC_COMPOSITE_ADDITIVE) { 322 pc->ops->apply = PCApply_Composite_Additive; 323 pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 324 } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 325 pc->ops->apply = PCApply_Composite_Multiplicative; 326 pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative; 327 } else if (type == PC_COMPOSITE_SPECIAL) { 328 pc->ops->apply = PCApply_Composite_Special; 329 pc->ops->applytranspose = NULL; 330 } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Unknown composite preconditioner type"); 331 jac->type = type; 332 PetscFunctionReturn(PETSC_SUCCESS); 333 } 334 335 static PetscErrorCode PCCompositeGetType_Composite(PC pc, PCCompositeType *type) 336 { 337 PC_Composite *jac = (PC_Composite *)pc->data; 338 339 PetscFunctionBegin; 340 *type = jac->type; 341 PetscFunctionReturn(PETSC_SUCCESS); 342 } 343 344 static PetscErrorCode PCCompositeAddPC_Composite(PC pc, PC subpc) 345 { 346 PC_Composite *jac; 347 PC_CompositeLink next, ilink; 348 PetscInt cnt = 0; 349 const char *prefix; 350 char newprefix[20]; 351 352 PetscFunctionBegin; 353 PetscCall(PetscNew(&ilink)); 354 ilink->next = NULL; 355 ilink->pc = subpc; 356 357 jac = (PC_Composite *)pc->data; 358 next = jac->head; 359 if (!next) { 360 jac->head = ilink; 361 ilink->previous = NULL; 362 } else { 363 cnt++; 364 while (next->next) { 365 next = next->next; 366 cnt++; 367 } 368 next->next = ilink; 369 ilink->previous = next; 370 } 371 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 372 PetscCall(PCSetOptionsPrefix(subpc, prefix)); 373 PetscCall(PetscSNPrintf(newprefix, 20, "sub_%d_", (int)cnt)); 374 PetscCall(PCAppendOptionsPrefix(subpc, newprefix)); 375 PetscCall(PetscObjectReference((PetscObject)subpc)); 376 PetscFunctionReturn(PETSC_SUCCESS); 377 } 378 379 static PetscErrorCode PCCompositeAddPCType_Composite(PC pc, PCType type) 380 { 381 PC subpc; 382 383 PetscFunctionBegin; 384 PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), &subpc)); 385 PetscCall(PetscObjectIncrementTabLevel((PetscObject)subpc, (PetscObject)pc, 1)); 386 PetscCall(PCCompositeAddPC_Composite(pc, subpc)); 387 /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 388 PetscCall(PCSetType(subpc, type)); 389 PetscCall(PCDestroy(&subpc)); 390 PetscFunctionReturn(PETSC_SUCCESS); 391 } 392 393 static PetscErrorCode PCCompositeGetNumberPC_Composite(PC pc, PetscInt *n) 394 { 395 PC_Composite *jac; 396 PC_CompositeLink next; 397 398 PetscFunctionBegin; 399 jac = (PC_Composite *)pc->data; 400 next = jac->head; 401 *n = 0; 402 while (next) { 403 next = next->next; 404 (*n)++; 405 } 406 PetscFunctionReturn(PETSC_SUCCESS); 407 } 408 409 static PetscErrorCode PCCompositeGetPC_Composite(PC pc, PetscInt n, PC *subpc) 410 { 411 PC_Composite *jac; 412 PC_CompositeLink next; 413 PetscInt i; 414 415 PetscFunctionBegin; 416 jac = (PC_Composite *)pc->data; 417 next = jac->head; 418 for (i = 0; i < n; i++) { 419 PetscCheck(next->next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Not enough PCs in composite preconditioner"); 420 next = next->next; 421 } 422 *subpc = next->pc; 423 PetscFunctionReturn(PETSC_SUCCESS); 424 } 425 426 /*@ 427 PCCompositeSetType - Sets the type of composite preconditioner. 428 429 Logically Collective 430 431 Input Parameters: 432 + pc - the preconditioner context 433 - type - `PC_COMPOSITE_ADDITIVE` (default), `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL` 434 435 Options Database Key: 436 . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 437 438 Level: advanced 439 440 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`, 441 `PCCompositeGetType()` 442 @*/ 443 PetscErrorCode PCCompositeSetType(PC pc, PCCompositeType type) 444 { 445 PetscFunctionBegin; 446 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 447 PetscValidLogicalCollectiveEnum(pc, type, 2); 448 PetscTryMethod(pc, "PCCompositeSetType_C", (PC, PCCompositeType), (pc, type)); 449 PetscFunctionReturn(PETSC_SUCCESS); 450 } 451 452 /*@ 453 PCCompositeGetType - Gets the type of composite preconditioner. 454 455 Logically Collective 456 457 Input Parameter: 458 . pc - the preconditioner context 459 460 Output Parameter: 461 . type - `PC_COMPOSITE_ADDITIVE` (default), `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL` 462 463 Level: advanced 464 465 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`, 466 `PCCompositeSetType()` 467 @*/ 468 PetscErrorCode PCCompositeGetType(PC pc, PCCompositeType *type) 469 { 470 PetscFunctionBegin; 471 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 472 PetscUseMethod(pc, "PCCompositeGetType_C", (PC, PCCompositeType *), (pc, type)); 473 PetscFunctionReturn(PETSC_SUCCESS); 474 } 475 476 /*@ 477 PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner, `PC_COMPOSITE_SPECIAL`, 478 for $\alpha I + R + S$ 479 480 Logically Collective 481 482 Input Parameters: 483 + pc - the preconditioner context 484 - alpha - scale on identity 485 486 Level: developer 487 488 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`, 489 `PCCompositeSetType()`, `PCCompositeGetType()` 490 @*/ 491 PetscErrorCode PCCompositeSpecialSetAlpha(PC pc, PetscScalar alpha) 492 { 493 PetscFunctionBegin; 494 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 495 PetscValidLogicalCollectiveScalar(pc, alpha, 2); 496 PetscTryMethod(pc, "PCCompositeSpecialSetAlpha_C", (PC, PetscScalar), (pc, alpha)); 497 PetscFunctionReturn(PETSC_SUCCESS); 498 } 499 500 PetscErrorCode PCCompositeSpecialSetAlphaMat(PC pc, Mat alpha_mat) 501 { 502 PetscFunctionBegin; 503 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 504 PetscTryMethod(pc, "PCCompositeSpecialSetAlphaMat_C", (PC, Mat), (pc, alpha_mat)); 505 PetscFunctionReturn(PETSC_SUCCESS); 506 } 507 508 /*@C 509 PCCompositeAddPCType - Adds another `PC` of the given type to the composite `PC`. 510 511 Collective 512 513 Input Parameters: 514 + pc - the preconditioner context 515 - type - the type of the new preconditioner 516 517 Level: intermediate 518 519 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeAddPC()`, `PCCompositeGetNumberPC()` 520 @*/ 521 PetscErrorCode PCCompositeAddPCType(PC pc, PCType type) 522 { 523 PetscFunctionBegin; 524 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 525 PetscTryMethod(pc, "PCCompositeAddPCType_C", (PC, PCType), (pc, type)); 526 PetscFunctionReturn(PETSC_SUCCESS); 527 } 528 529 /*@ 530 PCCompositeAddPC - Adds another `PC` to the composite `PC`. 531 532 Collective 533 534 Input Parameters: 535 + pc - the preconditioner context 536 - subpc - the new preconditioner 537 538 Level: intermediate 539 540 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeAddPCType()`, `PCCompositeGetNumberPC()` 541 @*/ 542 PetscErrorCode PCCompositeAddPC(PC pc, PC subpc) 543 { 544 PetscFunctionBegin; 545 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 546 PetscValidHeaderSpecific(subpc, PC_CLASSID, 2); 547 PetscTryMethod(pc, "PCCompositeAddPC_C", (PC, PC), (pc, subpc)); 548 PetscFunctionReturn(PETSC_SUCCESS); 549 } 550 551 /*@ 552 PCCompositeGetNumberPC - Gets the number of `PC` objects in the composite `PC`. 553 554 Not Collective 555 556 Input Parameter: 557 . pc - the preconditioner context 558 559 Output Parameter: 560 . num - the number of sub pcs 561 562 Level: developer 563 564 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeGetPC()`, `PCCompositeAddPC()`, `PCCompositeAddPCType()` 565 @*/ 566 PetscErrorCode PCCompositeGetNumberPC(PC pc, PetscInt *num) 567 { 568 PetscFunctionBegin; 569 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 570 PetscAssertPointer(num, 2); 571 PetscUseMethod(pc, "PCCompositeGetNumberPC_C", (PC, PetscInt *), (pc, num)); 572 PetscFunctionReturn(PETSC_SUCCESS); 573 } 574 575 /*@ 576 PCCompositeGetPC - Gets one of the `PC` objects in the composite `PC`. 577 578 Not Collective 579 580 Input Parameters: 581 + pc - the preconditioner context 582 - n - the number of the pc requested 583 584 Output Parameter: 585 . subpc - the PC requested 586 587 Level: intermediate 588 589 Note: 590 To use a different operator to construct one of the inner preconditioners first call `PCCompositeGetPC()`, then 591 call `PCSetOperators()` on that `PC`. 592 593 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeAddPCType()`, `PCCompositeGetNumberPC()`, `PCSetOperators()` 594 @*/ 595 PetscErrorCode PCCompositeGetPC(PC pc, PetscInt n, PC *subpc) 596 { 597 PetscFunctionBegin; 598 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 599 PetscAssertPointer(subpc, 3); 600 PetscUseMethod(pc, "PCCompositeGetPC_C", (PC, PetscInt, PC *), (pc, n, subpc)); 601 PetscFunctionReturn(PETSC_SUCCESS); 602 } 603 604 /*MC 605 PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 606 607 Options Database Keys: 608 + -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 609 . -pc_use_amat - activates `PCSetUseAmat()` 610 - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose 611 612 Level: intermediate 613 614 Notes: 615 To use a Krylov method inside the composite preconditioner, set the `PCType` of one or more 616 inner `PC`s to be `PCKSP`. Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 617 the incorrect answer) unless you use `KSPFGMRES` as the outer Krylov method 618 619 To use a different operator to construct one of the inner preconditioners first call `PCCompositeGetPC()`, then 620 call `PCSetOperators()` on that `PC`. 621 622 .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 623 `PCSHELL`, `PCKSP`, `PCCompositeSetType()`, `PCCompositeSpecialSetAlpha()`, `PCCompositeAddPCType()`, 624 `PCCompositeGetPC()`, `PCSetUseAmat()`, `PCCompositeAddPC()`, `PCCompositeGetNumberPC()` 625 M*/ 626 627 PETSC_EXTERN PetscErrorCode PCCreate_Composite(PC pc) 628 { 629 PC_Composite *jac; 630 631 PetscFunctionBegin; 632 PetscCall(PetscNew(&jac)); 633 634 pc->ops->apply = PCApply_Composite_Additive; 635 pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 636 pc->ops->setup = PCSetUp_Composite; 637 pc->ops->setuponblocks = PCSetUpOnBlocks_Composite; 638 pc->ops->reset = PCReset_Composite; 639 pc->ops->destroy = PCDestroy_Composite; 640 pc->ops->setfromoptions = PCSetFromOptions_Composite; 641 pc->ops->view = PCView_Composite; 642 pc->ops->applyrichardson = NULL; 643 644 pc->data = (void *)jac; 645 jac->type = PC_COMPOSITE_ADDITIVE; 646 jac->work1 = NULL; 647 jac->work2 = NULL; 648 jac->head = NULL; 649 jac->alpha_mat = NULL; 650 651 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSetType_C", PCCompositeSetType_Composite)); 652 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetType_C", PCCompositeGetType_Composite)); 653 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPCType_C", PCCompositeAddPCType_Composite)); 654 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPC_C", PCCompositeAddPC_Composite)); 655 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetNumberPC_C", PCCompositeGetNumberPC_Composite)); 656 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetPC_C", PCCompositeGetPC_Composite)); 657 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlpha_C", PCCompositeSpecialSetAlpha_Composite)); 658 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlphaMat_C", PCCompositeSpecialSetAlphaMat_Composite)); 659 PetscFunctionReturn(PETSC_SUCCESS); 660 } 661