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