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