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