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