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