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