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],*types[] = {"multiplicative","additive","special"}; 150 PetscTruth flg; 151 152 PetscFunctionBegin; 153 ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 154 ierr = PetscOptionsEList("-pc_composite_type","Type of composition","PCCompositeSetType",types,3,types[0],&indx,&flg);CHKERRQ(ierr); 155 if (flg) { 156 ierr = PCCompositeSetType(pc,(PCCompositeType)indx);CHKERRQ(ierr); 157 } 158 ierr = PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);CHKERRQ(ierr); 159 if (flg) { 160 ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr); 161 } 162 ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 163 if (flg) { 164 for (i=0; i<nmax; i++) { 165 ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 166 } 167 } 168 ierr = PetscOptionsTail();CHKERRQ(ierr); 169 170 next = jac->head; 171 while (next) { 172 ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 173 next = next->next; 174 } 175 PetscFunctionReturn(0); 176 } 177 178 #undef __FUNCT__ 179 #define __FUNCT__ "PCView_Composite" 180 static int PCView_Composite(PC pc,PetscViewer viewer) 181 { 182 PC_Composite *jac = (PC_Composite*)pc->data; 183 int ierr; 184 PC_CompositeLink next = jac->head; 185 PetscTruth isascii; 186 187 PetscFunctionBegin; 188 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 189 if (isascii) { 190 ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 191 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 192 } else { 193 SETERRQ1(1,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name); 194 } 195 if (isascii) { 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 (isascii) { 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 int 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 int 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(1,"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 int PCCompositeAddPC_Composite(PC pc,PCType type) 246 { 247 PC_Composite *jac; 248 PC_CompositeLink next,link; 249 int ierr,cnt = 0; 250 char *prefix,newprefix[8]; 251 252 PetscFunctionBegin; 253 ierr = PetscNew(struct _PC_CompositeLink,&link);CHKERRQ(ierr); 254 link->next = 0; 255 ierr = PCCreate(pc->comm,&link->pc);CHKERRQ(ierr); 256 257 jac = (PC_Composite*)pc->data; 258 next = jac->head; 259 if (!next) { 260 jac->head = link; 261 } else { 262 cnt++; 263 while (next->next) { 264 next = next->next; 265 cnt++; 266 } 267 next->next = link; 268 } 269 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 270 ierr = PCSetOptionsPrefix(link->pc,prefix);CHKERRQ(ierr); 271 sprintf(newprefix,"sub_%d_",cnt); 272 ierr = PCAppendOptionsPrefix(link->pc,newprefix);CHKERRQ(ierr); 273 /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 274 ierr = PCSetType(link->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 int PCCompositeGetPC_Composite(PC pc,int n,PC *subpc) 284 { 285 PC_Composite *jac; 286 PC_CompositeLink next; 287 int 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(1,"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 int 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 int PCCompositeSetType(PC pc,PCCompositeType type) 337 { 338 int ierr,(*f)(PC,PCCompositeType); 339 340 PetscFunctionBegin; 341 PetscValidHeaderSpecific(pc,PC_COOKIE); 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 int PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 366 { 367 int ierr,(*f)(PC,PetscScalar); 368 369 PetscFunctionBegin; 370 PetscValidHeaderSpecific(pc,PC_COOKIE); 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 int PCCompositeAddPC(PC pc,PCType type) 394 { 395 int ierr,(*f)(PC,PCType); 396 397 PetscFunctionBegin; 398 PetscValidHeaderSpecific(pc,PC_COOKIE); 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 int PCCompositeGetPC(PC pc,int n,PC *subpc) 427 { 428 int ierr,(*f)(PC,int,PC *); 429 430 PetscFunctionBegin; 431 PetscValidHeaderSpecific(pc,PC_COOKIE); 432 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 433 if (f) { 434 ierr = (*f)(pc,n,subpc);CHKERRQ(ierr); 435 } else { 436 SETERRQ(1,"Cannot get pc, not composite type"); 437 } 438 PetscFunctionReturn(0); 439 } 440 441 #undef __FUNCT__ 442 #define __FUNCT__ "PCCompositeSetUseTrue" 443 /*@ 444 PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than 445 the matrix used to define the preconditioner) is used to compute 446 the residual when the multiplicative scheme is used. 447 448 Collective on PC 449 450 Input Parameters: 451 . pc - the preconditioner context 452 453 Options Database Key: 454 . -pc_composite_true - Activates PCCompositeSetUseTrue() 455 456 Note: 457 For the common case in which the preconditioning and linear 458 system matrices are identical, this routine is unnecessary. 459 460 Level: Developer 461 462 .keywords: PC, composite preconditioner, set, true, flag 463 464 .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue() 465 @*/ 466 int PCCompositeSetUseTrue(PC pc) 467 { 468 int ierr,(*f)(PC); 469 470 PetscFunctionBegin; 471 PetscValidHeaderSpecific(pc,PC_COOKIE); 472 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr); 473 if (f) { 474 ierr = (*f)(pc);CHKERRQ(ierr); 475 } 476 PetscFunctionReturn(0); 477 } 478 479 /* -------------------------------------------------------------------------------------------*/ 480 481 /*MC 482 PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 483 484 Options Database Keys: 485 . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 486 . -pc_composite_true - Activates PCCompositeSetUseTrue() 487 488 Level: intermediate 489 490 Concepts: composing solvers 491 492 Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 493 inner PCs to be PCKSP. 494 Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 495 the incorrect answer) unless you use KSPFGMRES as the other Krylov method 496 497 498 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 499 PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 500 PCCompositeGetPC(), PCCompositeSetUseTrue() 501 502 M*/ 503 504 EXTERN_C_BEGIN 505 #undef __FUNCT__ 506 #define __FUNCT__ "PCCreate_Composite" 507 int PCCreate_Composite(PC pc) 508 { 509 int ierr; 510 PC_Composite *jac; 511 512 PetscFunctionBegin; 513 ierr = PetscNew(PC_Composite,&jac);CHKERRQ(ierr); 514 PetscLogObjectMemory(pc,sizeof(PC_Composite)); 515 pc->ops->apply = PCApply_Composite_Additive; 516 pc->ops->setup = PCSetUp_Composite; 517 pc->ops->destroy = PCDestroy_Composite; 518 pc->ops->setfromoptions = PCSetFromOptions_Composite; 519 pc->ops->view = PCView_Composite; 520 pc->ops->applyrichardson = 0; 521 522 pc->data = (void*)jac; 523 jac->type = PC_COMPOSITE_ADDITIVE; 524 jac->work1 = 0; 525 jac->work2 = 0; 526 jac->head = 0; 527 528 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite", 529 PCCompositeSetType_Composite);CHKERRQ(ierr); 530 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite", 531 PCCompositeAddPC_Composite);CHKERRQ(ierr); 532 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite", 533 PCCompositeGetPC_Composite);CHKERRQ(ierr); 534 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite", 535 PCCompositeSetUseTrue_Composite);CHKERRQ(ierr); 536 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite", 537 PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr); 538 539 PetscFunctionReturn(0); 540 } 541 EXTERN_C_END 542 543