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