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