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