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