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