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 PetscTruth 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; 130 131 PetscFunctionBegin; 132 while (next) { 133 ierr = PCDestroy(next->pc);CHKERRQ(ierr); 134 next = next->next; 135 } 136 137 if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);} 138 if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);} 139 ierr = PetscFree(jac);CHKERRQ(ierr); 140 PetscFunctionReturn(0); 141 } 142 143 #undef __FUNCT__ 144 #define __FUNCT__ "PCSetFromOptions_Composite" 145 static PetscErrorCode PCSetFromOptions_Composite(PC pc) 146 { 147 PC_Composite *jac = (PC_Composite*)pc->data; 148 PetscErrorCode ierr; 149 PetscInt nmax = 8,i; 150 PC_CompositeLink next; 151 char *pcs[8]; 152 PetscTruth flg; 153 154 PetscFunctionBegin; 155 ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 156 ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 157 if (flg) { 158 ierr = PCCompositeSetType(pc,jac->type);CHKERRQ(ierr); 159 } 160 flg = PETSC_FALSE; 161 ierr = PetscOptionsTruth("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 162 if (flg) { 163 ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr); 164 } 165 ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 166 if (flg) { 167 for (i=0; i<nmax; i++) { 168 ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 169 } 170 } 171 ierr = PetscOptionsTail();CHKERRQ(ierr); 172 173 next = jac->head; 174 while (next) { 175 ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 176 next = next->next; 177 } 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "PCView_Composite" 183 static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer) 184 { 185 PC_Composite *jac = (PC_Composite*)pc->data; 186 PetscErrorCode ierr; 187 PC_CompositeLink next = jac->head; 188 PetscTruth iascii; 189 190 PetscFunctionBegin; 191 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 192 if (iascii) { 193 ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr); 194 ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 195 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 196 } else { 197 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name); 198 } 199 if (iascii) { 200 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 201 } 202 while (next) { 203 ierr = PCView(next->pc,viewer);CHKERRQ(ierr); 204 next = next->next; 205 } 206 if (iascii) { 207 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 208 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 209 } 210 PetscFunctionReturn(0); 211 } 212 213 /* ------------------------------------------------------------------------------*/ 214 215 EXTERN_C_BEGIN 216 #undef __FUNCT__ 217 #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite" 218 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 219 { 220 PC_Composite *jac = (PC_Composite*)pc->data; 221 PetscFunctionBegin; 222 jac->alpha = alpha; 223 PetscFunctionReturn(0); 224 } 225 EXTERN_C_END 226 227 EXTERN_C_BEGIN 228 #undef __FUNCT__ 229 #define __FUNCT__ "PCCompositeSetType_Composite" 230 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType_Composite(PC pc,PCCompositeType type) 231 { 232 PetscFunctionBegin; 233 if (type == PC_COMPOSITE_ADDITIVE) { 234 pc->ops->apply = PCApply_Composite_Additive; 235 } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 236 pc->ops->apply = PCApply_Composite_Multiplicative; 237 } else if (type == PC_COMPOSITE_SPECIAL) { 238 pc->ops->apply = PCApply_Composite_Special; 239 } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type"); 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 const char *prefix; 254 char newprefix[8]; 255 256 PetscFunctionBegin; 257 ierr = PetscNewLog(pc,struct _PC_CompositeLink,&ilink);CHKERRQ(ierr); 258 ilink->next = 0; 259 ierr = PCCreate(((PetscObject)pc)->comm,&ilink->pc);CHKERRQ(ierr); 260 261 jac = (PC_Composite*)pc->data; 262 next = jac->head; 263 if (!next) { 264 jac->head = ilink; 265 ilink->previous = PETSC_NULL; 266 } else { 267 cnt++; 268 while (next->next) { 269 next = next->next; 270 cnt++; 271 } 272 next->next = ilink; 273 ilink->previous = next; 274 } 275 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 276 ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr); 277 sprintf(newprefix,"sub_%d_",(int)cnt); 278 ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr); 279 /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 280 ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr); 281 PetscFunctionReturn(0); 282 } 283 EXTERN_C_END 284 285 EXTERN_C_BEGIN 286 #undef __FUNCT__ 287 #define __FUNCT__ "PCCompositeGetPC_Composite" 288 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc) 289 { 290 PC_Composite *jac; 291 PC_CompositeLink next; 292 PetscInt i; 293 294 PetscFunctionBegin; 295 jac = (PC_Composite*)pc->data; 296 next = jac->head; 297 for (i=0; i<n; i++) { 298 if (!next->next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner"); 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 /*@ 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_CLASSID,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 /*@ 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_CLASSID,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_CLASSID,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 /*@ 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_CLASSID,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 SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type"); 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 PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC pc) 469 { 470 PetscErrorCode ierr,(*f)(PC); 471 472 PetscFunctionBegin; 473 PetscValidHeaderSpecific(pc,PC_CLASSID,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 - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose 490 491 Level: intermediate 492 493 Concepts: composing solvers 494 495 Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 496 inner PCs to be PCKSP. 497 Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 498 the incorrect answer) unless you use KSPFGMRES as the outter Krylov method 499 500 501 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 502 PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 503 PCCompositeGetPC(), PCCompositeSetUseTrue() 504 505 M*/ 506 507 EXTERN_C_BEGIN 508 #undef __FUNCT__ 509 #define __FUNCT__ "PCCreate_Composite" 510 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Composite(PC pc) 511 { 512 PetscErrorCode ierr; 513 PC_Composite *jac; 514 515 PetscFunctionBegin; 516 ierr = PetscNewLog(pc,PC_Composite,&jac);CHKERRQ(ierr); 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