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