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