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 DM dm; 172 173 PetscFunctionBegin; 174 if (!jac->work1) { 175 ierr = MatCreateVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr); 176 } 177 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 178 while (next) { 179 ierr = PCSetDM(next->pc,dm);CHKERRQ(ierr); 180 ierr = PCSetOperators(next->pc,pc->mat,pc->pmat);CHKERRQ(ierr); 181 next = next->next; 182 } 183 PetscFunctionReturn(0); 184 } 185 186 #undef __FUNCT__ 187 #define __FUNCT__ "PCReset_Composite" 188 static PetscErrorCode PCReset_Composite(PC pc) 189 { 190 PC_Composite *jac = (PC_Composite*)pc->data; 191 PetscErrorCode ierr; 192 PC_CompositeLink next = jac->head; 193 194 PetscFunctionBegin; 195 while (next) { 196 ierr = PCReset(next->pc);CHKERRQ(ierr); 197 next = next->next; 198 } 199 ierr = VecDestroy(&jac->work1);CHKERRQ(ierr); 200 ierr = VecDestroy(&jac->work2);CHKERRQ(ierr); 201 PetscFunctionReturn(0); 202 } 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "PCDestroy_Composite" 206 static PetscErrorCode PCDestroy_Composite(PC pc) 207 { 208 PC_Composite *jac = (PC_Composite*)pc->data; 209 PetscErrorCode ierr; 210 PC_CompositeLink next = jac->head,next_tmp; 211 212 PetscFunctionBegin; 213 ierr = PCReset_Composite(pc);CHKERRQ(ierr); 214 while (next) { 215 ierr = PCDestroy(&next->pc);CHKERRQ(ierr); 216 next_tmp = next; 217 next = next->next; 218 ierr = PetscFree(next_tmp);CHKERRQ(ierr); 219 } 220 ierr = PetscFree(pc->data);CHKERRQ(ierr); 221 PetscFunctionReturn(0); 222 } 223 224 #undef __FUNCT__ 225 #define __FUNCT__ "PCSetFromOptions_Composite" 226 static PetscErrorCode PCSetFromOptions_Composite(PetscOptions *PetscOptionsObject,PC pc) 227 { 228 PC_Composite *jac = (PC_Composite*)pc->data; 229 PetscErrorCode ierr; 230 PetscInt nmax = 8,i; 231 PC_CompositeLink next; 232 char *pcs[8]; 233 PetscBool flg; 234 235 PetscFunctionBegin; 236 ierr = PetscOptionsHead(PetscOptionsObject,"Composite preconditioner options");CHKERRQ(ierr); 237 ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 238 if (flg) { 239 ierr = PCCompositeSetType(pc,jac->type);CHKERRQ(ierr); 240 } 241 ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 242 if (flg) { 243 for (i=0; i<nmax; i++) { 244 ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 245 ierr = PetscFree(pcs[i]);CHKERRQ(ierr); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */ 246 } 247 } 248 ierr = PetscOptionsTail();CHKERRQ(ierr); 249 250 next = jac->head; 251 while (next) { 252 ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 253 next = next->next; 254 } 255 PetscFunctionReturn(0); 256 } 257 258 #undef __FUNCT__ 259 #define __FUNCT__ "PCView_Composite" 260 static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer) 261 { 262 PC_Composite *jac = (PC_Composite*)pc->data; 263 PetscErrorCode ierr; 264 PC_CompositeLink next = jac->head; 265 PetscBool iascii; 266 267 PetscFunctionBegin; 268 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 269 if (iascii) { 270 ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr); 271 ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 272 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 273 } 274 if (iascii) { 275 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 276 } 277 while (next) { 278 ierr = PCView(next->pc,viewer);CHKERRQ(ierr); 279 next = next->next; 280 } 281 if (iascii) { 282 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 283 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 284 } 285 PetscFunctionReturn(0); 286 } 287 288 /* ------------------------------------------------------------------------------*/ 289 290 #undef __FUNCT__ 291 #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite" 292 static PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 293 { 294 PC_Composite *jac = (PC_Composite*)pc->data; 295 296 PetscFunctionBegin; 297 jac->alpha = alpha; 298 PetscFunctionReturn(0); 299 } 300 301 #undef __FUNCT__ 302 #define __FUNCT__ "PCCompositeSetType_Composite" 303 static PetscErrorCode PCCompositeSetType_Composite(PC pc,PCCompositeType type) 304 { 305 PC_Composite *jac = (PC_Composite*)pc->data; 306 307 PetscFunctionBegin; 308 if (type == PC_COMPOSITE_ADDITIVE) { 309 pc->ops->apply = PCApply_Composite_Additive; 310 pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 311 } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 312 pc->ops->apply = PCApply_Composite_Multiplicative; 313 pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative; 314 } else if (type == PC_COMPOSITE_SPECIAL) { 315 pc->ops->apply = PCApply_Composite_Special; 316 pc->ops->applytranspose = NULL; 317 } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type"); 318 jac->type = type; 319 PetscFunctionReturn(0); 320 } 321 322 #undef __FUNCT__ 323 #define __FUNCT__ "PCCompositeGetType_Composite" 324 static PetscErrorCode PCCompositeGetType_Composite(PC pc,PCCompositeType *type) 325 { 326 PC_Composite *jac = (PC_Composite*)pc->data; 327 328 PetscFunctionBegin; 329 *type = jac->type; 330 PetscFunctionReturn(0); 331 } 332 333 #undef __FUNCT__ 334 #define __FUNCT__ "PCCompositeAddPC_Composite" 335 static PetscErrorCode PCCompositeAddPC_Composite(PC pc,PCType type) 336 { 337 PC_Composite *jac; 338 PC_CompositeLink next,ilink; 339 PetscErrorCode ierr; 340 PetscInt cnt = 0; 341 const char *prefix; 342 char newprefix[8]; 343 344 PetscFunctionBegin; 345 ierr = PetscNewLog(pc,&ilink);CHKERRQ(ierr); 346 ilink->next = 0; 347 ierr = PCCreate(PetscObjectComm((PetscObject)pc),&ilink->pc);CHKERRQ(ierr); 348 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->pc);CHKERRQ(ierr); 349 350 jac = (PC_Composite*)pc->data; 351 next = jac->head; 352 if (!next) { 353 jac->head = ilink; 354 ilink->previous = NULL; 355 } else { 356 cnt++; 357 while (next->next) { 358 next = next->next; 359 cnt++; 360 } 361 next->next = ilink; 362 ilink->previous = next; 363 } 364 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 365 ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr); 366 sprintf(newprefix,"sub_%d_",(int)cnt); 367 ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr); 368 /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 369 ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr); 370 PetscFunctionReturn(0); 371 } 372 373 #undef __FUNCT__ 374 #define __FUNCT__ "PCCompositeGetPC_Composite" 375 static PetscErrorCode PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc) 376 { 377 PC_Composite *jac; 378 PC_CompositeLink next; 379 PetscInt i; 380 381 PetscFunctionBegin; 382 jac = (PC_Composite*)pc->data; 383 next = jac->head; 384 for (i=0; i<n; i++) { 385 if (!next->next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner"); 386 next = next->next; 387 } 388 *subpc = next->pc; 389 PetscFunctionReturn(0); 390 } 391 392 /* -------------------------------------------------------------------------------- */ 393 #undef __FUNCT__ 394 #define __FUNCT__ "PCCompositeSetType" 395 /*@ 396 PCCompositeSetType - Sets the type of composite preconditioner. 397 398 Logically Collective on PC 399 400 Input Parameters: 401 + pc - the preconditioner context 402 - type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 403 404 Options Database Key: 405 . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 406 407 Level: Developer 408 409 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 410 @*/ 411 PetscErrorCode PCCompositeSetType(PC pc,PCCompositeType type) 412 { 413 PetscErrorCode ierr; 414 415 PetscFunctionBegin; 416 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 417 PetscValidLogicalCollectiveEnum(pc,type,2); 418 ierr = PetscTryMethod(pc,"PCCompositeSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 419 PetscFunctionReturn(0); 420 } 421 422 #undef __FUNCT__ 423 #define __FUNCT__ "PCCompositeGetType" 424 /*@ 425 PCCompositeGetType - Gets the type of composite preconditioner. 426 427 Logically Collective on PC 428 429 Input Parameter: 430 . pc - the preconditioner context 431 432 Output Parameter: 433 . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 434 435 Options Database Key: 436 . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 437 438 Level: Developer 439 440 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 441 @*/ 442 PetscErrorCode PCCompositeGetType(PC pc,PCCompositeType *type) 443 { 444 PetscErrorCode ierr; 445 446 PetscFunctionBegin; 447 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 448 ierr = PetscUseMethod(pc,"PCCompositeGetType_C",(PC,PCCompositeType*),(pc,type));CHKERRQ(ierr); 449 PetscFunctionReturn(0); 450 } 451 452 #undef __FUNCT__ 453 #define __FUNCT__ "PCCompositeSpecialSetAlpha" 454 /*@ 455 PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner 456 for alphaI + R + S 457 458 Logically Collective on PC 459 460 Input Parameter: 461 + pc - the preconditioner context 462 - alpha - scale on identity 463 464 Level: Developer 465 466 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 467 @*/ 468 PetscErrorCode PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 469 { 470 PetscErrorCode ierr; 471 472 PetscFunctionBegin; 473 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 474 PetscValidLogicalCollectiveScalar(pc,alpha,2); 475 ierr = PetscTryMethod(pc,"PCCompositeSpecialSetAlpha_C",(PC,PetscScalar),(pc,alpha));CHKERRQ(ierr); 476 PetscFunctionReturn(0); 477 } 478 479 #undef __FUNCT__ 480 #define __FUNCT__ "PCCompositeAddPC" 481 /*@C 482 PCCompositeAddPC - Adds another PC to the composite PC. 483 484 Collective on PC 485 486 Input Parameters: 487 + pc - the preconditioner context 488 - type - the type of the new preconditioner 489 490 Level: Developer 491 492 .keywords: PC, composite preconditioner, add 493 @*/ 494 PetscErrorCode PCCompositeAddPC(PC pc,PCType type) 495 { 496 PetscErrorCode ierr; 497 498 PetscFunctionBegin; 499 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 500 ierr = PetscTryMethod(pc,"PCCompositeAddPC_C",(PC,PCType),(pc,type));CHKERRQ(ierr); 501 PetscFunctionReturn(0); 502 } 503 504 #undef __FUNCT__ 505 #define __FUNCT__ "PCCompositeGetPC" 506 /*@ 507 PCCompositeGetPC - Gets one of the PC objects in the composite PC. 508 509 Not Collective 510 511 Input Parameter: 512 + pc - the preconditioner context 513 - n - the number of the pc requested 514 515 Output Parameters: 516 . subpc - the PC requested 517 518 Level: Developer 519 520 .keywords: PC, get, composite preconditioner, sub preconditioner 521 522 .seealso: PCCompositeAddPC() 523 @*/ 524 PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *subpc) 525 { 526 PetscErrorCode ierr; 527 528 PetscFunctionBegin; 529 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 530 PetscValidPointer(subpc,3); 531 ierr = PetscUseMethod(pc,"PCCompositeGetPC_C",(PC,PetscInt,PC*),(pc,n,subpc));CHKERRQ(ierr); 532 PetscFunctionReturn(0); 533 } 534 535 /* -------------------------------------------------------------------------------------------*/ 536 537 /*MC 538 PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 539 540 Options Database Keys: 541 + -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 542 . -pc_use_amat - Activates PCSetUseAmat() 543 - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose 544 545 Level: intermediate 546 547 Concepts: composing solvers 548 549 Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 550 inner PCs to be PCKSP. 551 Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 552 the incorrect answer) unless you use KSPFGMRES as the outer Krylov method 553 554 555 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 556 PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 557 PCCompositeGetPC(), PCSetUseAmat() 558 559 M*/ 560 561 #undef __FUNCT__ 562 #define __FUNCT__ "PCCreate_Composite" 563 PETSC_EXTERN PetscErrorCode PCCreate_Composite(PC pc) 564 { 565 PetscErrorCode ierr; 566 PC_Composite *jac; 567 568 PetscFunctionBegin; 569 ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 570 571 pc->ops->apply = PCApply_Composite_Additive; 572 pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 573 pc->ops->setup = PCSetUp_Composite; 574 pc->ops->reset = PCReset_Composite; 575 pc->ops->destroy = PCDestroy_Composite; 576 pc->ops->setfromoptions = PCSetFromOptions_Composite; 577 pc->ops->view = PCView_Composite; 578 pc->ops->applyrichardson = 0; 579 580 pc->data = (void*)jac; 581 jac->type = PC_COMPOSITE_ADDITIVE; 582 jac->work1 = 0; 583 jac->work2 = 0; 584 jac->head = 0; 585 586 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSetType_C",PCCompositeSetType_Composite);CHKERRQ(ierr); 587 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetType_C",PCCompositeGetType_Composite);CHKERRQ(ierr); 588 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCCompositeAddPC_C",PCCompositeAddPC_Composite);CHKERRQ(ierr); 589 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetPC_C",PCCompositeGetPC_Composite);CHKERRQ(ierr); 590 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr); 591 PetscFunctionReturn(0); 592 } 593 594