1 2 /* 3 Defines a preconditioner that can consist of a collection of PCs 4 */ 5 #include <private/pcimpl.h> /*I "petscpc.h" I*/ 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 PetscBool use_true_matrix; 22 } PC_Composite; 23 24 #undef __FUNCT__ 25 #define __FUNCT__ "PCApply_Composite_Multiplicative" 26 static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y) 27 { 28 PetscErrorCode ierr; 29 PC_Composite *jac = (PC_Composite*)pc->data; 30 PC_CompositeLink next = jac->head; 31 Mat mat = pc->pmat; 32 33 PetscFunctionBegin; 34 if (!next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 35 if (next->next && !jac->work2) { /* allocate second work vector */ 36 ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr); 37 } 38 if (jac->use_true_matrix) mat = pc->mat; 39 ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 40 while (next->next) { 41 next = next->next; 42 ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr); 43 ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr); 44 ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 45 ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 46 ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 47 } 48 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 49 while (next->previous) { 50 next = next->previous; 51 ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr); 52 ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr); 53 ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 54 ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 55 ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 56 } 57 } 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "PCApplyTranspose_Composite_Multiplicative" 63 static PetscErrorCode PCApplyTranspose_Composite_Multiplicative(PC pc,Vec x,Vec y) 64 { 65 PetscErrorCode ierr; 66 PC_Composite *jac = (PC_Composite*)pc->data; 67 PC_CompositeLink next = jac->head; 68 Mat mat = pc->pmat; 69 70 PetscFunctionBegin; 71 if (!next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 72 if (next->next && !jac->work2) { /* allocate second work vector */ 73 ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr); 74 } 75 if (jac->use_true_matrix) mat = pc->mat; 76 /* locate last PC */ 77 while (next->next) { 78 next = next->next; 79 } 80 ierr = PCApplyTranspose(next->pc,x,y);CHKERRQ(ierr); 81 while (next->previous) { 82 next = next->previous; 83 ierr = MatMultTranspose(mat,y,jac->work1);CHKERRQ(ierr); 84 ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr); 85 ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 86 ierr = PCApplyTranspose(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 87 ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 88 } 89 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 90 next = jac->head; 91 while (next->next) { 92 next = next->next; 93 ierr = MatMultTranspose(mat,y,jac->work1);CHKERRQ(ierr); 94 ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr); 95 ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 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 #undef __FUNCT__ 109 #define __FUNCT__ "PCApply_Composite_Special" 110 static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y) 111 { 112 PetscErrorCode ierr; 113 PC_Composite *jac = (PC_Composite*)pc->data; 114 PC_CompositeLink next = jac->head; 115 116 PetscFunctionBegin; 117 if (!next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 118 if (!next->next || next->next->next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs"); 119 120 ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 121 ierr = PCApply(next->next->pc,jac->work1,y);CHKERRQ(ierr); 122 PetscFunctionReturn(0); 123 } 124 125 #undef __FUNCT__ 126 #define __FUNCT__ "PCApply_Composite_Additive" 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(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 135 ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 136 while (next->next) { 137 next = next->next; 138 ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 139 ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 140 ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 141 } 142 PetscFunctionReturn(0); 143 } 144 145 #undef __FUNCT__ 146 #define __FUNCT__ "PCApplyTranspose_Composite_Additive" 147 static PetscErrorCode PCApplyTranspose_Composite_Additive(PC pc,Vec x,Vec y) 148 { 149 PetscErrorCode ierr; 150 PC_Composite *jac = (PC_Composite*)pc->data; 151 PC_CompositeLink next = jac->head; 152 153 PetscFunctionBegin; 154 if (!next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 155 ierr = PCApplyTranspose(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 = PCApplyTranspose(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__ "PCSetUp_Composite" 167 static PetscErrorCode PCSetUp_Composite(PC pc) 168 { 169 PetscErrorCode ierr; 170 PC_Composite *jac = (PC_Composite*)pc->data; 171 PC_CompositeLink next = jac->head; 172 173 PetscFunctionBegin; 174 if (!jac->work1) { 175 ierr = MatGetVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr); 176 } 177 while (next) { 178 ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 179 next = next->next; 180 } 181 PetscFunctionReturn(0); 182 } 183 184 #undef __FUNCT__ 185 #define __FUNCT__ "PCReset_Composite" 186 static PetscErrorCode PCReset_Composite(PC pc) 187 { 188 PC_Composite *jac = (PC_Composite*)pc->data; 189 PetscErrorCode ierr; 190 PC_CompositeLink next = jac->head; 191 192 PetscFunctionBegin; 193 while (next) { 194 ierr = PCReset(next->pc);CHKERRQ(ierr); 195 next = next->next; 196 } 197 ierr = VecDestroy(&jac->work1);CHKERRQ(ierr); 198 ierr = VecDestroy(&jac->work2);CHKERRQ(ierr); 199 PetscFunctionReturn(0); 200 } 201 202 #undef __FUNCT__ 203 #define __FUNCT__ "PCDestroy_Composite" 204 static PetscErrorCode PCDestroy_Composite(PC pc) 205 { 206 PC_Composite *jac = (PC_Composite*)pc->data; 207 PetscErrorCode ierr; 208 PC_CompositeLink next = jac->head,next_tmp; 209 210 PetscFunctionBegin; 211 ierr = PCReset_Composite(pc);CHKERRQ(ierr); 212 while (next) { 213 ierr = PCDestroy(&next->pc);CHKERRQ(ierr); 214 next_tmp = next; 215 next = next->next; 216 ierr = PetscFree(next_tmp);CHKERRQ(ierr); 217 } 218 ierr = PetscFree(pc->data);CHKERRQ(ierr); 219 PetscFunctionReturn(0); 220 } 221 222 #undef __FUNCT__ 223 #define __FUNCT__ "PCSetFromOptions_Composite" 224 static PetscErrorCode PCSetFromOptions_Composite(PC pc) 225 { 226 PC_Composite *jac = (PC_Composite*)pc->data; 227 PetscErrorCode ierr; 228 PetscInt nmax = 8,i; 229 PC_CompositeLink next; 230 char *pcs[8]; 231 PetscBool flg; 232 233 PetscFunctionBegin; 234 ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 235 ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 236 if (flg) { 237 ierr = PCCompositeSetType(pc,jac->type);CHKERRQ(ierr); 238 } 239 flg = PETSC_FALSE; 240 ierr = PetscOptionsBool("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 241 if (flg) { 242 ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr); 243 } 244 ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 245 if (flg) { 246 for (i=0; i<nmax; i++) { 247 ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 248 ierr = PetscFree(pcs[i]);CHKERRQ(ierr); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */ 249 } 250 } 251 ierr = PetscOptionsTail();CHKERRQ(ierr); 252 253 next = jac->head; 254 while (next) { 255 ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 256 next = next->next; 257 } 258 PetscFunctionReturn(0); 259 } 260 261 #undef __FUNCT__ 262 #define __FUNCT__ "PCView_Composite" 263 static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer) 264 { 265 PC_Composite *jac = (PC_Composite*)pc->data; 266 PetscErrorCode ierr; 267 PC_CompositeLink next = jac->head; 268 PetscBool iascii; 269 270 PetscFunctionBegin; 271 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 272 if (iascii) { 273 ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr); 274 ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 275 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 276 } else { 277 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name); 278 } 279 if (iascii) { 280 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 281 } 282 while (next) { 283 ierr = PCView(next->pc,viewer);CHKERRQ(ierr); 284 next = next->next; 285 } 286 if (iascii) { 287 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 288 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 289 } 290 PetscFunctionReturn(0); 291 } 292 293 /* ------------------------------------------------------------------------------*/ 294 295 EXTERN_C_BEGIN 296 #undef __FUNCT__ 297 #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite" 298 PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 299 { 300 PC_Composite *jac = (PC_Composite*)pc->data; 301 PetscFunctionBegin; 302 jac->alpha = alpha; 303 PetscFunctionReturn(0); 304 } 305 EXTERN_C_END 306 307 EXTERN_C_BEGIN 308 #undef __FUNCT__ 309 #define __FUNCT__ "PCCompositeSetType_Composite" 310 PetscErrorCode PCCompositeSetType_Composite(PC pc,PCCompositeType type) 311 { 312 PetscFunctionBegin; 313 if (type == PC_COMPOSITE_ADDITIVE) { 314 pc->ops->apply = PCApply_Composite_Additive; 315 pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 316 } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 317 pc->ops->apply = PCApply_Composite_Multiplicative; 318 pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative; 319 } else if (type == PC_COMPOSITE_SPECIAL) { 320 pc->ops->apply = PCApply_Composite_Special; 321 pc->ops->applytranspose = PETSC_NULL; 322 } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type"); 323 PetscFunctionReturn(0); 324 } 325 EXTERN_C_END 326 327 EXTERN_C_BEGIN 328 #undef __FUNCT__ 329 #define __FUNCT__ "PCCompositeAddPC_Composite" 330 PetscErrorCode PCCompositeAddPC_Composite(PC pc,PCType type) 331 { 332 PC_Composite *jac; 333 PC_CompositeLink next,ilink; 334 PetscErrorCode ierr; 335 PetscInt cnt = 0; 336 const char *prefix; 337 char newprefix[8]; 338 339 PetscFunctionBegin; 340 ierr = PetscNewLog(pc,struct _PC_CompositeLink,&ilink);CHKERRQ(ierr); 341 ilink->next = 0; 342 ierr = PCCreate(((PetscObject)pc)->comm,&ilink->pc);CHKERRQ(ierr); 343 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->pc);CHKERRQ(ierr); 344 345 jac = (PC_Composite*)pc->data; 346 next = jac->head; 347 if (!next) { 348 jac->head = ilink; 349 ilink->previous = PETSC_NULL; 350 } else { 351 cnt++; 352 while (next->next) { 353 next = next->next; 354 cnt++; 355 } 356 next->next = ilink; 357 ilink->previous = next; 358 } 359 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 360 ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr); 361 sprintf(newprefix,"sub_%d_",(int)cnt); 362 ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr); 363 /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 364 ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr); 365 PetscFunctionReturn(0); 366 } 367 EXTERN_C_END 368 369 EXTERN_C_BEGIN 370 #undef __FUNCT__ 371 #define __FUNCT__ "PCCompositeGetPC_Composite" 372 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(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner"); 383 next = next->next; 384 } 385 *subpc = next->pc; 386 PetscFunctionReturn(0); 387 } 388 EXTERN_C_END 389 390 EXTERN_C_BEGIN 391 #undef __FUNCT__ 392 #define __FUNCT__ "PCCompositeSetUseTrue_Composite" 393 PetscErrorCode PCCompositeSetUseTrue_Composite(PC pc) 394 { 395 PC_Composite *jac; 396 397 PetscFunctionBegin; 398 jac = (PC_Composite*)pc->data; 399 jac->use_true_matrix = PETSC_TRUE; 400 PetscFunctionReturn(0); 401 } 402 EXTERN_C_END 403 404 /* -------------------------------------------------------------------------------- */ 405 #undef __FUNCT__ 406 #define __FUNCT__ "PCCompositeSetType" 407 /*@ 408 PCCompositeSetType - Sets the type of composite preconditioner. 409 410 Logically Collective on PC 411 412 Input Parameter: 413 + pc - the preconditioner context 414 - type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 415 416 Options Database Key: 417 . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 418 419 Level: Developer 420 421 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 422 @*/ 423 PetscErrorCode PCCompositeSetType(PC pc,PCCompositeType type) 424 { 425 PetscErrorCode ierr; 426 427 PetscFunctionBegin; 428 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 429 PetscValidLogicalCollectiveEnum(pc,type,2); 430 ierr = PetscTryMethod(pc,"PCCompositeSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 431 PetscFunctionReturn(0); 432 } 433 434 #undef __FUNCT__ 435 #define __FUNCT__ "PCCompositeSpecialSetAlpha" 436 /*@ 437 PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner 438 for alphaI + R + S 439 440 Logically Collective on PC 441 442 Input Parameter: 443 + pc - the preconditioner context 444 - alpha - scale on identity 445 446 Level: Developer 447 448 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 449 @*/ 450 PetscErrorCode PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 451 { 452 PetscErrorCode ierr; 453 454 PetscFunctionBegin; 455 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 456 PetscValidLogicalCollectiveScalar(pc,alpha,2); 457 ierr = PetscTryMethod(pc,"PCCompositeSpecialSetAlpha_C",(PC,PetscScalar),(pc,alpha));CHKERRQ(ierr); 458 PetscFunctionReturn(0); 459 } 460 461 #undef __FUNCT__ 462 #define __FUNCT__ "PCCompositeAddPC" 463 /*@C 464 PCCompositeAddPC - Adds another PC to the composite PC. 465 466 Collective on PC 467 468 Input Parameters: 469 + pc - the preconditioner context 470 - type - the type of the new preconditioner 471 472 Level: Developer 473 474 .keywords: PC, composite preconditioner, add 475 @*/ 476 PetscErrorCode PCCompositeAddPC(PC pc,PCType type) 477 { 478 PetscErrorCode ierr; 479 480 PetscFunctionBegin; 481 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 482 ierr = PetscTryMethod(pc,"PCCompositeAddPC_C",(PC,PCType),(pc,type));CHKERRQ(ierr); 483 PetscFunctionReturn(0); 484 } 485 486 #undef __FUNCT__ 487 #define __FUNCT__ "PCCompositeGetPC" 488 /*@ 489 PCCompositeGetPC - Gets one of the PC objects in the composite PC. 490 491 Not Collective 492 493 Input Parameter: 494 + pc - the preconditioner context 495 - n - the number of the pc requested 496 497 Output Parameters: 498 . subpc - the PC requested 499 500 Level: Developer 501 502 .keywords: PC, get, composite preconditioner, sub preconditioner 503 504 .seealso: PCCompositeAddPC() 505 @*/ 506 PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *subpc) 507 { 508 PetscErrorCode ierr; 509 510 PetscFunctionBegin; 511 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 512 PetscValidPointer(subpc,3); 513 ierr = PetscUseMethod(pc,"PCCompositeGetPC_C",(PC,PetscInt,PC *),(pc,n,subpc));CHKERRQ(ierr); 514 PetscFunctionReturn(0); 515 } 516 517 #undef __FUNCT__ 518 #define __FUNCT__ "PCCompositeSetUseTrue" 519 /*@ 520 PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than 521 the matrix used to define the preconditioner) is used to compute 522 the residual when the multiplicative scheme is used. 523 524 Logically Collective on PC 525 526 Input Parameters: 527 . pc - the preconditioner context 528 529 Options Database Key: 530 . -pc_composite_true - Activates PCCompositeSetUseTrue() 531 532 Note: 533 For the common case in which the preconditioning and linear 534 system matrices are identical, this routine is unnecessary. 535 536 Level: Developer 537 538 .keywords: PC, composite preconditioner, set, true, flag 539 540 .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue() 541 @*/ 542 PetscErrorCode PCCompositeSetUseTrue(PC pc) 543 { 544 PetscErrorCode ierr; 545 546 PetscFunctionBegin; 547 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 548 ierr = PetscTryMethod(pc,"PCCompositeSetUseTrue_C",(PC),(pc));CHKERRQ(ierr); 549 PetscFunctionReturn(0); 550 } 551 552 /* -------------------------------------------------------------------------------------------*/ 553 554 /*MC 555 PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 556 557 Options Database Keys: 558 + -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 559 . -pc_composite_true - Activates PCCompositeSetUseTrue() 560 - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose 561 562 Level: intermediate 563 564 Concepts: composing solvers 565 566 Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 567 inner PCs to be PCKSP. 568 Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 569 the incorrect answer) unless you use KSPFGMRES as the outter Krylov method 570 571 572 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 573 PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 574 PCCompositeGetPC(), PCCompositeSetUseTrue() 575 576 M*/ 577 578 EXTERN_C_BEGIN 579 #undef __FUNCT__ 580 #define __FUNCT__ "PCCreate_Composite" 581 PetscErrorCode PCCreate_Composite(PC pc) 582 { 583 PetscErrorCode ierr; 584 PC_Composite *jac; 585 586 PetscFunctionBegin; 587 ierr = PetscNewLog(pc,PC_Composite,&jac);CHKERRQ(ierr); 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 = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite", 604 PCCompositeSetType_Composite);CHKERRQ(ierr); 605 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite", 606 PCCompositeAddPC_Composite);CHKERRQ(ierr); 607 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite", 608 PCCompositeGetPC_Composite);CHKERRQ(ierr); 609 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite", 610 PCCompositeSetUseTrue_Composite);CHKERRQ(ierr); 611 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite", 612 PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr); 613 614 PetscFunctionReturn(0); 615 } 616 EXTERN_C_END 617 618