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