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(PCGetFailedReason(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, i; 243 PC_CompositeLink next; 244 char *pcs[1024]; 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 nmax = (PetscInt)PETSC_STATIC_ARRAY_LENGTH(pcs); 252 PetscCall(PetscOptionsStringArray("-pc_composite_pcs", "List of composite solvers", "PCCompositeAddPCType", pcs, &nmax, &flg)); 253 if (flg) { 254 for (i = 0; i < nmax; i++) { 255 PetscCall(PCCompositeAddPCType(pc, pcs[i])); 256 PetscCall(PetscFree(pcs[i])); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */ 257 } 258 } 259 PetscOptionsHeadEnd(); 260 261 next = jac->head; 262 while (next) { 263 PetscCall(PCSetFromOptions(next->pc)); 264 next = next->next; 265 } 266 PetscFunctionReturn(PETSC_SUCCESS); 267 } 268 269 static PetscErrorCode PCView_Composite(PC pc, PetscViewer viewer) 270 { 271 PC_Composite *jac = (PC_Composite *)pc->data; 272 PC_CompositeLink next = jac->head; 273 PetscBool iascii; 274 275 PetscFunctionBegin; 276 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 277 if (iascii) { 278 PetscCall(PetscViewerASCIIPrintf(viewer, "Composite PC type - %s\n", PCCompositeTypes[jac->type])); 279 PetscCall(PetscViewerASCIIPrintf(viewer, "PCs on composite preconditioner follow\n")); 280 PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------\n")); 281 } 282 if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer)); 283 while (next) { 284 PetscCall(PCView(next->pc, viewer)); 285 next = next->next; 286 } 287 if (iascii) { 288 PetscCall(PetscViewerASCIIPopTab(viewer)); 289 PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------\n")); 290 } 291 PetscFunctionReturn(PETSC_SUCCESS); 292 } 293 294 static PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc, PetscScalar alpha) 295 { 296 PC_Composite *jac = (PC_Composite *)pc->data; 297 298 PetscFunctionBegin; 299 jac->alpha = alpha; 300 PetscFunctionReturn(PETSC_SUCCESS); 301 } 302 303 static PetscErrorCode PCCompositeSpecialSetAlphaMat_Composite(PC pc, Mat alpha_mat) 304 { 305 PC_Composite *jac = (PC_Composite *)pc->data; 306 307 PetscFunctionBegin; 308 if (alpha_mat) { 309 PetscValidHeaderSpecific(alpha_mat, MAT_CLASSID, 2); 310 PetscCall(PetscObjectReference((PetscObject)alpha_mat)); 311 } 312 PetscCall(MatDestroy(&jac->alpha_mat)); 313 jac->alpha_mat = alpha_mat; 314 PetscFunctionReturn(PETSC_SUCCESS); 315 } 316 317 static PetscErrorCode PCCompositeSetType_Composite(PC pc, PCCompositeType type) 318 { 319 PC_Composite *jac = (PC_Composite *)pc->data; 320 321 PetscFunctionBegin; 322 if (type == PC_COMPOSITE_ADDITIVE) { 323 pc->ops->apply = PCApply_Composite_Additive; 324 pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 325 } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 326 pc->ops->apply = PCApply_Composite_Multiplicative; 327 pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative; 328 } else if (type == PC_COMPOSITE_SPECIAL) { 329 pc->ops->apply = PCApply_Composite_Special; 330 pc->ops->applytranspose = NULL; 331 } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Unknown composite preconditioner type"); 332 jac->type = type; 333 PetscFunctionReturn(PETSC_SUCCESS); 334 } 335 336 static PetscErrorCode PCCompositeGetType_Composite(PC pc, PCCompositeType *type) 337 { 338 PC_Composite *jac = (PC_Composite *)pc->data; 339 340 PetscFunctionBegin; 341 *type = jac->type; 342 PetscFunctionReturn(PETSC_SUCCESS); 343 } 344 345 static PetscErrorCode PCCompositeAddPC_Composite(PC pc, PC subpc) 346 { 347 PC_Composite *jac; 348 PC_CompositeLink next, ilink; 349 PetscInt cnt = 0; 350 const char *prefix; 351 char newprefix[20]; 352 353 PetscFunctionBegin; 354 PetscCall(PetscNew(&ilink)); 355 ilink->next = NULL; 356 ilink->pc = subpc; 357 358 jac = (PC_Composite *)pc->data; 359 next = jac->head; 360 if (!next) { 361 jac->head = ilink; 362 ilink->previous = NULL; 363 } else { 364 cnt++; 365 while (next->next) { 366 next = next->next; 367 cnt++; 368 } 369 next->next = ilink; 370 ilink->previous = next; 371 } 372 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 373 PetscCall(PCSetOptionsPrefix(subpc, prefix)); 374 PetscCall(PetscSNPrintf(newprefix, 20, "sub_%" PetscInt_FMT "_", cnt)); 375 PetscCall(PCAppendOptionsPrefix(subpc, newprefix)); 376 PetscCall(PetscObjectReference((PetscObject)subpc)); 377 PetscFunctionReturn(PETSC_SUCCESS); 378 } 379 380 static PetscErrorCode PCCompositeAddPCType_Composite(PC pc, PCType type) 381 { 382 PC subpc; 383 384 PetscFunctionBegin; 385 PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), &subpc)); 386 PetscCall(PetscObjectIncrementTabLevel((PetscObject)subpc, (PetscObject)pc, 1)); 387 PetscCall(PCCompositeAddPC_Composite(pc, subpc)); 388 /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 389 PetscCall(PCSetType(subpc, type)); 390 PetscCall(PCDestroy(&subpc)); 391 PetscFunctionReturn(PETSC_SUCCESS); 392 } 393 394 static PetscErrorCode PCCompositeGetNumberPC_Composite(PC pc, PetscInt *n) 395 { 396 PC_Composite *jac; 397 PC_CompositeLink next; 398 399 PetscFunctionBegin; 400 jac = (PC_Composite *)pc->data; 401 next = jac->head; 402 *n = 0; 403 while (next) { 404 next = next->next; 405 (*n)++; 406 } 407 PetscFunctionReturn(PETSC_SUCCESS); 408 } 409 410 static PetscErrorCode PCCompositeGetPC_Composite(PC pc, PetscInt n, PC *subpc) 411 { 412 PC_Composite *jac; 413 PC_CompositeLink next; 414 PetscInt i; 415 416 PetscFunctionBegin; 417 jac = (PC_Composite *)pc->data; 418 next = jac->head; 419 for (i = 0; i < n; i++) { 420 PetscCheck(next->next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Not enough PCs in composite preconditioner"); 421 next = next->next; 422 } 423 *subpc = next->pc; 424 PetscFunctionReturn(PETSC_SUCCESS); 425 } 426 427 /*@ 428 PCCompositeSetType - Sets the type of composite preconditioner. 429 430 Logically Collective 431 432 Input Parameters: 433 + pc - the preconditioner context 434 - type - `PC_COMPOSITE_ADDITIVE` (default), `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL` 435 436 Options Database Key: 437 . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 438 439 Level: advanced 440 441 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`, 442 `PCCompositeGetType()` 443 @*/ 444 PetscErrorCode PCCompositeSetType(PC pc, PCCompositeType type) 445 { 446 PetscFunctionBegin; 447 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 448 PetscValidLogicalCollectiveEnum(pc, type, 2); 449 PetscTryMethod(pc, "PCCompositeSetType_C", (PC, PCCompositeType), (pc, type)); 450 PetscFunctionReturn(PETSC_SUCCESS); 451 } 452 453 /*@ 454 PCCompositeGetType - Gets the type of composite preconditioner. 455 456 Logically Collective 457 458 Input Parameter: 459 . pc - the preconditioner context 460 461 Output Parameter: 462 . type - `PC_COMPOSITE_ADDITIVE` (default), `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL` 463 464 Level: advanced 465 466 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`, 467 `PCCompositeSetType()` 468 @*/ 469 PetscErrorCode PCCompositeGetType(PC pc, PCCompositeType *type) 470 { 471 PetscFunctionBegin; 472 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 473 PetscUseMethod(pc, "PCCompositeGetType_C", (PC, PCCompositeType *), (pc, type)); 474 PetscFunctionReturn(PETSC_SUCCESS); 475 } 476 477 /*@ 478 PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner, `PC_COMPOSITE_SPECIAL`, 479 for $\alpha I + R + S$ 480 481 Logically Collective 482 483 Input Parameters: 484 + pc - the preconditioner context 485 - alpha - scale on identity 486 487 Level: developer 488 489 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`, 490 `PCCompositeSetType()`, `PCCompositeGetType()` 491 @*/ 492 PetscErrorCode PCCompositeSpecialSetAlpha(PC pc, PetscScalar alpha) 493 { 494 PetscFunctionBegin; 495 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 496 PetscValidLogicalCollectiveScalar(pc, alpha, 2); 497 PetscTryMethod(pc, "PCCompositeSpecialSetAlpha_C", (PC, PetscScalar), (pc, alpha)); 498 PetscFunctionReturn(PETSC_SUCCESS); 499 } 500 501 PetscErrorCode PCCompositeSpecialSetAlphaMat(PC pc, Mat alpha_mat) 502 { 503 PetscFunctionBegin; 504 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 505 PetscTryMethod(pc, "PCCompositeSpecialSetAlphaMat_C", (PC, Mat), (pc, alpha_mat)); 506 PetscFunctionReturn(PETSC_SUCCESS); 507 } 508 509 /*@C 510 PCCompositeAddPCType - Adds another `PC` of the given type to the composite `PC`. 511 512 Collective 513 514 Input Parameters: 515 + pc - the preconditioner context 516 - type - the type of the new preconditioner 517 518 Level: intermediate 519 520 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeAddPC()`, `PCCompositeGetNumberPC()` 521 @*/ 522 PetscErrorCode PCCompositeAddPCType(PC pc, PCType type) 523 { 524 PetscFunctionBegin; 525 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 526 PetscTryMethod(pc, "PCCompositeAddPCType_C", (PC, PCType), (pc, type)); 527 PetscFunctionReturn(PETSC_SUCCESS); 528 } 529 530 /*@ 531 PCCompositeAddPC - Adds another `PC` to the composite `PC`. 532 533 Collective 534 535 Input Parameters: 536 + pc - the preconditioner context 537 - subpc - the new preconditioner 538 539 Level: intermediate 540 541 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeAddPCType()`, `PCCompositeGetNumberPC()` 542 @*/ 543 PetscErrorCode PCCompositeAddPC(PC pc, PC subpc) 544 { 545 PetscFunctionBegin; 546 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 547 PetscValidHeaderSpecific(subpc, PC_CLASSID, 2); 548 PetscTryMethod(pc, "PCCompositeAddPC_C", (PC, PC), (pc, subpc)); 549 PetscFunctionReturn(PETSC_SUCCESS); 550 } 551 552 /*@ 553 PCCompositeGetNumberPC - Gets the number of `PC` objects in the composite `PC`. 554 555 Not Collective 556 557 Input Parameter: 558 . pc - the preconditioner context 559 560 Output Parameter: 561 . num - the number of sub pcs 562 563 Level: developer 564 565 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeGetPC()`, `PCCompositeAddPC()`, `PCCompositeAddPCType()` 566 @*/ 567 PetscErrorCode PCCompositeGetNumberPC(PC pc, PetscInt *num) 568 { 569 PetscFunctionBegin; 570 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 571 PetscAssertPointer(num, 2); 572 PetscUseMethod(pc, "PCCompositeGetNumberPC_C", (PC, PetscInt *), (pc, num)); 573 PetscFunctionReturn(PETSC_SUCCESS); 574 } 575 576 /*@ 577 PCCompositeGetPC - Gets one of the `PC` objects in the composite `PC`. 578 579 Not Collective 580 581 Input Parameters: 582 + pc - the preconditioner context 583 - n - the number of the pc requested 584 585 Output Parameter: 586 . subpc - the PC requested 587 588 Level: intermediate 589 590 Note: 591 To use a different operator to construct one of the inner preconditioners first call `PCCompositeGetPC()`, then 592 call `PCSetOperators()` on that `PC`. 593 594 .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeAddPCType()`, `PCCompositeGetNumberPC()`, `PCSetOperators()` 595 @*/ 596 PetscErrorCode PCCompositeGetPC(PC pc, PetscInt n, PC *subpc) 597 { 598 PetscFunctionBegin; 599 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 600 PetscAssertPointer(subpc, 3); 601 PetscUseMethod(pc, "PCCompositeGetPC_C", (PC, PetscInt, PC *), (pc, n, subpc)); 602 PetscFunctionReturn(PETSC_SUCCESS); 603 } 604 605 /*MC 606 PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 607 608 Options Database Keys: 609 + -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 610 . -pc_use_amat - activates `PCSetUseAmat()` 611 - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose 612 613 Level: intermediate 614 615 Notes: 616 To use a Krylov method inside the composite preconditioner, set the `PCType` of one or more 617 inner `PC`s to be `PCKSP`. Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 618 the incorrect answer) unless you use `KSPFGMRES` as the outer Krylov method 619 620 To use a different operator to construct one of the inner preconditioners first call `PCCompositeGetPC()`, then 621 call `PCSetOperators()` on that `PC`. 622 623 .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 624 `PCSHELL`, `PCKSP`, `PCCompositeSetType()`, `PCCompositeSpecialSetAlpha()`, `PCCompositeAddPCType()`, 625 `PCCompositeGetPC()`, `PCSetUseAmat()`, `PCCompositeAddPC()`, `PCCompositeGetNumberPC()` 626 M*/ 627 628 PETSC_EXTERN PetscErrorCode PCCreate_Composite(PC pc) 629 { 630 PC_Composite *jac; 631 632 PetscFunctionBegin; 633 PetscCall(PetscNew(&jac)); 634 635 pc->ops->apply = PCApply_Composite_Additive; 636 pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 637 pc->ops->setup = PCSetUp_Composite; 638 pc->ops->setuponblocks = PCSetUpOnBlocks_Composite; 639 pc->ops->reset = PCReset_Composite; 640 pc->ops->destroy = PCDestroy_Composite; 641 pc->ops->setfromoptions = PCSetFromOptions_Composite; 642 pc->ops->view = PCView_Composite; 643 pc->ops->applyrichardson = NULL; 644 645 pc->data = (void *)jac; 646 jac->type = PC_COMPOSITE_ADDITIVE; 647 jac->work1 = NULL; 648 jac->work2 = NULL; 649 jac->head = NULL; 650 jac->alpha_mat = NULL; 651 652 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSetType_C", PCCompositeSetType_Composite)); 653 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetType_C", PCCompositeGetType_Composite)); 654 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPCType_C", PCCompositeAddPCType_Composite)); 655 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPC_C", PCCompositeAddPC_Composite)); 656 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetNumberPC_C", PCCompositeGetNumberPC_Composite)); 657 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetPC_C", PCCompositeGetPC_Composite)); 658 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlpha_C", PCCompositeSpecialSetAlpha_Composite)); 659 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlphaMat_C", PCCompositeSpecialSetAlphaMat_Composite)); 660 PetscFunctionReturn(PETSC_SUCCESS); 661 } 662