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