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