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