1 /*$Id: composite.c,v 1.45 2001/08/07 03:03:39 balay Exp $*/ 2 /* 3 Defines a preconditioner that can consist of a collection of PCs 4 */ 5 #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 6 #include "petscksp.h" /*I "petscksp.h" I*/ 7 8 typedef struct _PC_CompositeLink *PC_CompositeLink; 9 struct _PC_CompositeLink { 10 PC pc; 11 PC_CompositeLink next; 12 }; 13 14 typedef struct { 15 PC_CompositeLink head; 16 PCCompositeType type; 17 Vec work1; 18 Vec work2; 19 PetscScalar alpha; 20 PetscTruth use_true_matrix; 21 } PC_Composite; 22 23 #undef __FUNCT__ 24 #define __FUNCT__ "PCApply_Composite_Multiplicative" 25 static int PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y) 26 { 27 int ierr; 28 PC_Composite *jac = (PC_Composite*)pc->data; 29 PC_CompositeLink next = jac->head; 30 PetscScalar one = 1.0,mone = -1.0; 31 Mat mat = pc->pmat; 32 33 PetscFunctionBegin; 34 if (!next) { 35 SETERRQ(1,"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,PC_LEFT);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(&mone,jac->work1,x,jac->work2);CHKERRQ(ierr); 46 ierr = PCApply(next->pc,jac->work2,jac->work1,PC_LEFT);CHKERRQ(ierr); 47 ierr = VecAXPY(&one,jac->work1,y);CHKERRQ(ierr); 48 } 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 int PCApply_Composite_Special(PC pc,Vec x,Vec y) 61 { 62 int ierr; 63 PC_Composite *jac = (PC_Composite*)pc->data; 64 PC_CompositeLink next = jac->head; 65 66 PetscFunctionBegin; 67 if (!next) { 68 SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 69 } 70 if (!next->next || next->next->next) { 71 SETERRQ(1,"Special composite preconditioners requires exactly two PCs"); 72 } 73 74 ierr = PCApply(next->pc,x,jac->work1,PC_LEFT);CHKERRQ(ierr); 75 ierr = PCApply(next->next->pc,jac->work1,y,PC_LEFT);CHKERRQ(ierr); 76 PetscFunctionReturn(0); 77 } 78 79 #undef __FUNCT__ 80 #define __FUNCT__ "PCApply_Composite_Additive" 81 static int PCApply_Composite_Additive(PC pc,Vec x,Vec y) 82 { 83 int ierr; 84 PC_Composite *jac = (PC_Composite*)pc->data; 85 PC_CompositeLink next = jac->head; 86 PetscScalar one = 1.0; 87 88 PetscFunctionBegin; 89 if (!next) { 90 SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 91 } 92 ierr = PCApply(next->pc,x,y,PC_LEFT);CHKERRQ(ierr); 93 while (next->next) { 94 next = next->next; 95 ierr = PCApply(next->pc,x,jac->work1,PC_LEFT);CHKERRQ(ierr); 96 ierr = VecAXPY(&one,jac->work1,y);CHKERRQ(ierr); 97 } 98 PetscFunctionReturn(0); 99 } 100 101 #undef __FUNCT__ 102 #define __FUNCT__ "PCSetUp_Composite" 103 static int PCSetUp_Composite(PC pc) 104 { 105 int ierr; 106 PC_Composite *jac = (PC_Composite*)pc->data; 107 PC_CompositeLink next = jac->head; 108 109 PetscFunctionBegin; 110 if (!jac->work1) { 111 ierr = VecDuplicate(pc->vec,&jac->work1);CHKERRQ(ierr); 112 } 113 while (next) { 114 ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 115 ierr = PCSetVector(next->pc,jac->work1);CHKERRQ(ierr); 116 next = next->next; 117 } 118 119 PetscFunctionReturn(0); 120 } 121 122 #undef __FUNCT__ 123 #define __FUNCT__ "PCDestroy_Composite" 124 static int PCDestroy_Composite(PC pc) 125 { 126 PC_Composite *jac = (PC_Composite*)pc->data; 127 int 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 int PCSetFromOptions_Composite(PC pc) 145 { 146 PC_Composite *jac = (PC_Composite*)pc->data; 147 int ierr,nmax = 8,i,indx; 148 PC_CompositeLink next; 149 char *pcs[8]; 150 const char *types[] = {"multiplicative","additive","special"}; 151 PetscTruth flg; 152 153 PetscFunctionBegin; 154 ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 155 ierr = PetscOptionsEList("-pc_composite_type","Type of composition","PCCompositeSetType",types,3,types[0],&indx,&flg);CHKERRQ(ierr); 156 if (flg) { 157 ierr = PCCompositeSetType(pc,(PCCompositeType)indx);CHKERRQ(ierr); 158 } 159 ierr = PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);CHKERRQ(ierr); 160 if (flg) { 161 ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr); 162 } 163 ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 164 if (flg) { 165 for (i=0; i<nmax; i++) { 166 ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 167 } 168 } 169 ierr = PetscOptionsTail();CHKERRQ(ierr); 170 171 next = jac->head; 172 while (next) { 173 ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 174 next = next->next; 175 } 176 PetscFunctionReturn(0); 177 } 178 179 #undef __FUNCT__ 180 #define __FUNCT__ "PCView_Composite" 181 static int PCView_Composite(PC pc,PetscViewer viewer) 182 { 183 PC_Composite *jac = (PC_Composite*)pc->data; 184 int ierr; 185 PC_CompositeLink next = jac->head; 186 PetscTruth isascii; 187 188 PetscFunctionBegin; 189 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 190 if (isascii) { 191 ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 192 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 193 } else { 194 SETERRQ1(1,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name); 195 } 196 if (isascii) { 197 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 198 } 199 while (next) { 200 ierr = PCView(next->pc,viewer);CHKERRQ(ierr); 201 next = next->next; 202 } 203 if (isascii) { 204 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 205 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 206 } 207 PetscFunctionReturn(0); 208 } 209 210 /* ------------------------------------------------------------------------------*/ 211 212 EXTERN_C_BEGIN 213 #undef __FUNCT__ 214 #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite" 215 int PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 216 { 217 PC_Composite *jac = (PC_Composite*)pc->data; 218 PetscFunctionBegin; 219 jac->alpha = alpha; 220 PetscFunctionReturn(0); 221 } 222 EXTERN_C_END 223 224 EXTERN_C_BEGIN 225 #undef __FUNCT__ 226 #define __FUNCT__ "PCCompositeSetType_Composite" 227 int PCCompositeSetType_Composite(PC pc,PCCompositeType type) 228 { 229 PetscFunctionBegin; 230 if (type == PC_COMPOSITE_ADDITIVE) { 231 pc->ops->apply = PCApply_Composite_Additive; 232 } else if (type == PC_COMPOSITE_MULTIPLICATIVE) { 233 pc->ops->apply = PCApply_Composite_Multiplicative; 234 } else if (type == PC_COMPOSITE_SPECIAL) { 235 pc->ops->apply = PCApply_Composite_Special; 236 } else { 237 SETERRQ(1,"Unkown composite preconditioner type"); 238 } 239 PetscFunctionReturn(0); 240 } 241 EXTERN_C_END 242 243 EXTERN_C_BEGIN 244 #undef __FUNCT__ 245 #define __FUNCT__ "PCCompositeAddPC_Composite" 246 int PCCompositeAddPC_Composite(PC pc,PCType type) 247 { 248 PC_Composite *jac; 249 PC_CompositeLink next,link; 250 int ierr,cnt = 0; 251 char *prefix,newprefix[8]; 252 253 PetscFunctionBegin; 254 ierr = PetscNew(struct _PC_CompositeLink,&link);CHKERRQ(ierr); 255 link->next = 0; 256 ierr = PCCreate(pc->comm,&link->pc);CHKERRQ(ierr); 257 258 jac = (PC_Composite*)pc->data; 259 next = jac->head; 260 if (!next) { 261 jac->head = link; 262 } else { 263 cnt++; 264 while (next->next) { 265 next = next->next; 266 cnt++; 267 } 268 next->next = link; 269 } 270 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 271 ierr = PCSetOptionsPrefix(link->pc,prefix);CHKERRQ(ierr); 272 sprintf(newprefix,"sub_%d_",cnt); 273 ierr = PCAppendOptionsPrefix(link->pc,newprefix);CHKERRQ(ierr); 274 /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 275 ierr = PCSetType(link->pc,type);CHKERRQ(ierr); 276 277 PetscFunctionReturn(0); 278 } 279 EXTERN_C_END 280 281 EXTERN_C_BEGIN 282 #undef __FUNCT__ 283 #define __FUNCT__ "PCCompositeGetPC_Composite" 284 int PCCompositeGetPC_Composite(PC pc,int n,PC *subpc) 285 { 286 PC_Composite *jac; 287 PC_CompositeLink next; 288 int 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(1,"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 int 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 /*@C 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 int PCCompositeSetType(PC pc,PCCompositeType type) 338 { 339 int 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 /*@C 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 int PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 367 { 368 int 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 int PCCompositeAddPC(PC pc,PCType type) 395 { 396 int 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 /*@C 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 int PCCompositeGetPC(PC pc,int n,PC *subpc) 428 { 429 int ierr,(*f)(PC,int,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(1,"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 int PCCompositeSetUseTrue(PC pc) 469 { 470 int 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 490 Level: intermediate 491 492 Concepts: composing solvers 493 494 Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 495 inner PCs to be PCKSP. 496 Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 497 the incorrect answer) unless you use KSPFGMRES as the other Krylov method 498 499 500 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 501 PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 502 PCCompositeGetPC(), PCCompositeSetUseTrue() 503 504 M*/ 505 506 EXTERN_C_BEGIN 507 #undef __FUNCT__ 508 #define __FUNCT__ "PCCreate_Composite" 509 int PCCreate_Composite(PC pc) 510 { 511 int ierr; 512 PC_Composite *jac; 513 514 PetscFunctionBegin; 515 ierr = PetscNew(PC_Composite,&jac);CHKERRQ(ierr); 516 PetscLogObjectMemory(pc,sizeof(PC_Composite)); 517 pc->ops->apply = PCApply_Composite_Additive; 518 pc->ops->setup = PCSetUp_Composite; 519 pc->ops->destroy = PCDestroy_Composite; 520 pc->ops->setfromoptions = PCSetFromOptions_Composite; 521 pc->ops->view = PCView_Composite; 522 pc->ops->applyrichardson = 0; 523 524 pc->data = (void*)jac; 525 jac->type = PC_COMPOSITE_ADDITIVE; 526 jac->work1 = 0; 527 jac->work2 = 0; 528 jac->head = 0; 529 530 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite", 531 PCCompositeSetType_Composite);CHKERRQ(ierr); 532 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite", 533 PCCompositeAddPC_Composite);CHKERRQ(ierr); 534 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite", 535 PCCompositeGetPC_Composite);CHKERRQ(ierr); 536 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite", 537 PCCompositeSetUseTrue_Composite);CHKERRQ(ierr); 538 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite", 539 PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr); 540 541 PetscFunctionReturn(0); 542 } 543 EXTERN_C_END 544 545