1 #define PETSCKSP_DLL 2 3 /* 4 Defines the multigrid preconditioner interface. 5 */ 6 #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 7 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "PCMGMCycle_Private" 11 PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG **mglevels,PetscTruth *converged) 12 { 13 PC_MG *mg = *mglevels,*mgc; 14 PetscErrorCode ierr; 15 PetscInt cycles = (PetscInt) mg->cycles; 16 17 PetscFunctionBegin; 18 if (converged) *converged = PETSC_FALSE; 19 20 if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 21 ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr); /* pre-smooth */ 22 if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 23 if (mg->level) { /* not the coarsest grid */ 24 if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);} 25 ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr); 26 if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);} 27 28 /* if on finest level and have convergence criteria set */ 29 if (mg->level == mg->levels-1 && mg->ttol) { 30 PetscReal rnorm; 31 ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr); 32 if (rnorm <= mg->ttol) { 33 *converged = PETSC_TRUE; 34 if (rnorm < mg->abstol) { 35 ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr); 36 } else { 37 ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr); 38 } 39 PetscFunctionReturn(0); 40 } 41 } 42 43 mgc = *(mglevels - 1); 44 if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 45 ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr); 46 if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 47 ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 48 while (cycles--) { 49 ierr = PCMGMCycle_Private(pc,mglevels-1,converged);CHKERRQ(ierr); 50 } 51 if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 52 ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr); 53 if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 54 if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 55 ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr); /* post smooth */ 56 if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 57 } 58 PetscFunctionReturn(0); 59 } 60 61 /* 62 PCMGCreate_Private - Creates a PC_MG structure for use with the 63 multigrid code. Level 0 is the coarsest. (But the 64 finest level is stored first in the array). 65 66 */ 67 #undef __FUNCT__ 68 #define __FUNCT__ "PCMGCreate_Private" 69 static PetscErrorCode PCMGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,PC_MG ***result) 70 { 71 PC_MG **mg; 72 PetscErrorCode ierr; 73 PetscInt i; 74 PetscMPIInt size; 75 const char *prefix; 76 PC ipc; 77 78 PetscFunctionBegin; 79 ierr = PetscMalloc(levels*sizeof(PC_MG*),&mg);CHKERRQ(ierr); 80 ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 81 82 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 83 84 for (i=0; i<levels; i++) { 85 ierr = PetscNewLog(pc,PC_MG,&mg[i]);CHKERRQ(ierr); 86 mg[i]->level = i; 87 mg[i]->levels = levels; 88 mg[i]->cycles = PC_MG_CYCLE_V; 89 mg[i]->galerkin = PETSC_FALSE; 90 mg[i]->galerkinused = PETSC_FALSE; 91 mg[i]->default_smoothu = 1; 92 mg[i]->default_smoothd = 1; 93 94 if (comms) comm = comms[i]; 95 ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr); 96 ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr); 97 ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr); 98 99 /* do special stuff for coarse grid */ 100 if (!i && levels > 1) { 101 ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 102 103 /* coarse solve is (redundant) LU by default */ 104 ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 105 ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr); 106 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 107 if (size > 1) { 108 ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 109 } else { 110 ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 111 } 112 113 } else { 114 char tprefix[128]; 115 sprintf(tprefix,"mg_levels_%d_",(int)i); 116 ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr); 117 } 118 ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr); 119 mg[i]->smoothu = mg[i]->smoothd; 120 mg[i]->rtol = 0.0; 121 mg[i]->abstol = 0.0; 122 mg[i]->dtol = 0.0; 123 mg[i]->ttol = 0.0; 124 mg[i]->eventsmoothsetup = 0; 125 mg[i]->eventsmoothsolve = 0; 126 mg[i]->eventresidual = 0; 127 mg[i]->eventinterprestrict = 0; 128 mg[i]->cyclesperpcapply = 1; 129 } 130 *result = mg; 131 PetscFunctionReturn(0); 132 } 133 134 #undef __FUNCT__ 135 #define __FUNCT__ "PCDestroy_MG" 136 static PetscErrorCode PCDestroy_MG(PC pc) 137 { 138 PC_MG **mg = (PC_MG**)pc->data; 139 PetscErrorCode ierr; 140 PetscInt i,n; 141 142 PetscFunctionBegin; 143 if (!mg) PetscFunctionReturn(0); 144 n = mg[0]->levels; 145 for (i=0; i<n-1; i++) { 146 if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);} 147 if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);} 148 if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);} 149 if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);} 150 if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);} 151 } 152 153 for (i=0; i<n; i++) { 154 if (mg[i]->smoothd != mg[i]->smoothu) { 155 ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr); 156 } 157 ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr); 158 ierr = PetscFree(mg[i]);CHKERRQ(ierr); 159 } 160 ierr = PetscFree(mg);CHKERRQ(ierr); 161 PetscFunctionReturn(0); 162 } 163 164 165 166 EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**); 167 EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG**); 168 EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**); 169 170 /* 171 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 172 or full cycle of multigrid. 173 174 Note: 175 A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 176 */ 177 #undef __FUNCT__ 178 #define __FUNCT__ "PCApply_MG" 179 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 180 { 181 PC_MG **mg = (PC_MG**)pc->data; 182 PetscErrorCode ierr; 183 PetscInt levels = mg[0]->levels,i; 184 185 PetscFunctionBegin; 186 mg[levels-1]->b = b; 187 mg[levels-1]->x = x; 188 if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) { 189 Vec tvec; 190 ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr); 191 ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr); 192 ierr = VecDestroy(tvec);CHKERRQ(ierr); 193 } 194 if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 195 ierr = VecSet(x,0.0);CHKERRQ(ierr); 196 for (i=0; i<mg[0]->cyclesperpcapply; i++) { 197 ierr = PCMGMCycle_Private(pc,mg+levels-1,PETSC_NULL);CHKERRQ(ierr); 198 } 199 } 200 else if (mg[0]->am == PC_MG_ADDITIVE) { 201 ierr = PCMGACycle_Private(mg);CHKERRQ(ierr); 202 } 203 else if (mg[0]->am == PC_MG_KASKADE) { 204 ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr); 205 } 206 else { 207 ierr = PCMGFCycle_Private(pc,mg);CHKERRQ(ierr); 208 } 209 PetscFunctionReturn(0); 210 } 211 212 #undef __FUNCT__ 213 #define __FUNCT__ "PCApplyRichardson_MG" 214 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 215 { 216 PC_MG **mg = (PC_MG**)pc->data; 217 PetscErrorCode ierr; 218 PetscInt levels = mg[0]->levels; 219 PetscTruth converged = PETSC_FALSE; 220 221 PetscFunctionBegin; 222 mg[levels-1]->b = b; 223 mg[levels-1]->x = x; 224 225 mg[levels-1]->rtol = rtol; 226 mg[levels-1]->abstol = abstol; 227 mg[levels-1]->dtol = dtol; 228 if (rtol) { 229 /* compute initial residual norm for relative convergence test */ 230 PetscReal rnorm; 231 ierr = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr); 232 ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 233 mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol); 234 } else if (abstol) { 235 mg[levels-1]->ttol = abstol; 236 } else { 237 mg[levels-1]->ttol = 0.0; 238 } 239 240 while (its-- && !converged) { 241 ierr = PCMGMCycle_Private(pc,mg+levels-1,&converged);CHKERRQ(ierr); 242 } 243 PetscFunctionReturn(0); 244 } 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "PCSetFromOptions_MG" 248 PetscErrorCode PCSetFromOptions_MG(PC pc) 249 { 250 PetscErrorCode ierr; 251 PetscInt m,levels = 1,cycles; 252 PetscTruth flg; 253 PC_MG **mg = (PC_MG**)pc->data; 254 PCMGType mgtype = PC_MG_ADDITIVE; 255 PCMGCycleType mgctype; 256 257 PetscFunctionBegin; 258 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 259 if (!pc->data) { 260 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 261 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 262 mg = (PC_MG**)pc->data; 263 } 264 mgctype = (PCMGCycleType) mg[0]->cycles; 265 ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 266 if (flg) { 267 ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 268 }; 269 ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr); 270 if (flg) { 271 ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 272 } 273 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 274 if (flg) { 275 ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 276 } 277 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 278 if (flg) { 279 ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 280 } 281 ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 282 if (flg) { 283 ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 284 } 285 if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 286 ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 287 if (flg) { 288 ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 289 } 290 } 291 ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr); 292 if (flg) { 293 PetscInt i; 294 char eventname[128]; 295 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 296 levels = mg[0]->levels; 297 for (i=0; i<levels; i++) { 298 sprintf(eventname,"MGSetup Level %d",(int)i); 299 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsetup);CHKERRQ(ierr); 300 sprintf(eventname,"MGSmooth Level %d",(int)i); 301 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsolve);CHKERRQ(ierr); 302 if (i) { 303 sprintf(eventname,"MGResid Level %d",(int)i); 304 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventresidual);CHKERRQ(ierr); 305 sprintf(eventname,"MGInterp Level %d",(int)i); 306 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventinterprestrict);CHKERRQ(ierr); 307 } 308 } 309 } 310 ierr = PetscOptionsTail();CHKERRQ(ierr); 311 PetscFunctionReturn(0); 312 } 313 314 const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 315 const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 316 317 #undef __FUNCT__ 318 #define __FUNCT__ "PCView_MG" 319 static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 320 { 321 PC_MG **mg = (PC_MG**)pc->data; 322 PetscErrorCode ierr; 323 PetscInt levels = mg[0]->levels,i; 324 PetscTruth iascii; 325 326 PetscFunctionBegin; 327 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 328 if (iascii) { 329 ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr); 330 if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 331 ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg[0]->cyclesperpcapply);CHKERRQ(ierr); 332 } 333 if (mg[0]->galerkin) { 334 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 335 } 336 for (i=0; i<levels; i++) { 337 if (!i) { 338 ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D presmooths=%D postsmooths=%D -----\n",i,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr); 339 } else { 340 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothd);CHKERRQ(ierr); 341 } 342 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 343 ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr); 344 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 345 if (i && mg[i]->smoothd == mg[i]->smoothu) { 346 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 347 } else if (i){ 348 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothu);CHKERRQ(ierr); 349 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 350 ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr); 351 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 352 } 353 } 354 } else { 355 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 356 } 357 PetscFunctionReturn(0); 358 } 359 360 /* 361 Calls setup for the KSP on each level 362 */ 363 #undef __FUNCT__ 364 #define __FUNCT__ "PCSetUp_MG" 365 static PetscErrorCode PCSetUp_MG(PC pc) 366 { 367 PC_MG **mg = (PC_MG**)pc->data; 368 PetscErrorCode ierr; 369 PetscInt i,n = mg[0]->levels; 370 PC cpc,mpc; 371 PetscTruth preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset; 372 PetscViewerASCIIMonitor ascii; 373 PetscViewer viewer = PETSC_NULL; 374 MPI_Comm comm; 375 Mat dA,dB; 376 MatStructure uflag; 377 Vec tvec; 378 379 PetscFunctionBegin; 380 381 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 382 /* so use those from global PC */ 383 /* Is this what we always want? What if user wants to keep old one? */ 384 ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 385 ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 386 ierr = KSPGetPC(mg[n-1]->smoothd,&mpc);CHKERRQ(ierr); 387 if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2))) { 388 ierr = PetscInfo(pc,"Using outer operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr); 389 ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 390 } 391 392 if (mg[0]->galerkin) { 393 Mat B; 394 mg[0]->galerkinused = PETSC_TRUE; 395 /* currently only handle case where mat and pmat are the same on coarser levels */ 396 ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 397 if (!pc->setupcalled) { 398 for (i=n-2; i>-1; i--) { 399 ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 400 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 401 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 402 dB = B; 403 } 404 ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr); 405 } else { 406 for (i=n-2; i>-1; i--) { 407 ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 408 ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 409 ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 410 dB = B; 411 } 412 } 413 } 414 415 if (!pc->setupcalled) { 416 ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr); 417 418 for (i=0; i<n; i++) { 419 if (monitor) { 420 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr); 421 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 422 ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 423 } 424 ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr); 425 } 426 for (i=1; i<n; i++) { 427 if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) { 428 if (monitor) { 429 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr); 430 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 431 ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 432 } 433 ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr); 434 } 435 } 436 for (i=1; i<n; i++) { 437 if (!mg[i]->residual) { 438 Mat mat; 439 ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 440 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 441 } 442 if (mg[i]->restrct && !mg[i]->interpolate) { 443 ierr = PCMGSetInterpolation(pc,i,mg[i]->restrct);CHKERRQ(ierr); 444 } 445 if (!mg[i]->restrct && mg[i]->interpolate) { 446 ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr); 447 } 448 #if defined(PETSC_USE_DEBUG) 449 if (!mg[i]->restrct || !mg[i]->interpolate) { 450 SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 451 } 452 #endif 453 } 454 for (i=0; i<n-1; i++) { 455 if (!mg[i]->b) { 456 Vec *vec; 457 ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 458 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 459 ierr = PetscFree(vec);CHKERRQ(ierr); 460 } 461 if (!mg[i]->r && i) { 462 ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 463 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 464 ierr = VecDestroy(tvec);CHKERRQ(ierr); 465 } 466 if (!mg[i]->x) { 467 ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 468 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 469 ierr = VecDestroy(tvec);CHKERRQ(ierr); 470 } 471 } 472 if (n != 1 && !mg[n-1]->r) { 473 /* PCMGSetR() on the finest level if user did not supply it */ 474 Vec *vec; 475 ierr = KSPGetVecs(mg[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 476 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 477 ierr = PetscFree(vec);CHKERRQ(ierr); 478 } 479 } 480 481 482 for (i=1; i<n; i++) { 483 if (mg[i]->smoothu == mg[i]->smoothd) { 484 /* if doing only down then initial guess is zero */ 485 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 486 } 487 if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 488 ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 489 if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 490 } 491 for (i=1; i<n; i++) { 492 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 493 Mat downmat,downpmat; 494 MatStructure matflag; 495 PetscTruth opsset; 496 497 /* check if operators have been set for up, if not use down operators to set them */ 498 ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 499 if (!opsset) { 500 ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 501 ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 502 } 503 504 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 505 if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 506 ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 507 if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 508 } 509 } 510 511 /* 512 If coarse solver is not direct method then DO NOT USE preonly 513 */ 514 ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 515 if (preonly) { 516 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 517 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 518 ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 519 if (!lu && !redundant && !cholesky) { 520 ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 521 } 522 } 523 524 if (!pc->setupcalled) { 525 if (monitor) { 526 ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 527 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 528 ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 529 } 530 ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr); 531 } 532 533 if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 534 ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr); 535 if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 536 537 /* 538 Dump the interpolation/restriction matrices plus the 539 Jacobian/stiffness on each level. This allows Matlab users to 540 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 541 542 Only support one or the other at the same time. 543 */ 544 #if defined(PETSC_USE_SOCKET_VIEWER) 545 ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr); 546 if (dump) { 547 viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 548 } 549 #endif 550 ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr); 551 if (dump) { 552 viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 553 } 554 555 if (viewer) { 556 for (i=1; i<n; i++) { 557 ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr); 558 } 559 for (i=0; i<n; i++) { 560 ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 561 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 562 } 563 } 564 PetscFunctionReturn(0); 565 } 566 567 /* -------------------------------------------------------------------------------------*/ 568 569 #undef __FUNCT__ 570 #define __FUNCT__ "PCMGSetLevels" 571 /*@C 572 PCMGSetLevels - Sets the number of levels to use with MG. 573 Must be called before any other MG routine. 574 575 Collective on PC 576 577 Input Parameters: 578 + pc - the preconditioner context 579 . levels - the number of levels 580 - comms - optional communicators for each level; this is to allow solving the coarser problems 581 on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 582 583 Level: intermediate 584 585 Notes: 586 If the number of levels is one then the multigrid uses the -mg_levels prefix 587 for setting the level options rather than the -mg_coarse prefix. 588 589 .keywords: MG, set, levels, multigrid 590 591 .seealso: PCMGSetType(), PCMGGetLevels() 592 @*/ 593 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 594 { 595 PetscErrorCode ierr; 596 PC_MG **mg=0; 597 598 PetscFunctionBegin; 599 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 600 601 if (pc->data) { 602 SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\ 603 make sure that you call PCMGSetLevels() before KSPSetFromOptions()"); 604 } 605 ierr = PCMGCreate_Private(((PetscObject)pc)->comm,levels,pc,comms,&mg);CHKERRQ(ierr); 606 mg[0]->am = PC_MG_MULTIPLICATIVE; 607 pc->data = (void*)mg; 608 pc->ops->applyrichardson = PCApplyRichardson_MG; 609 PetscFunctionReturn(0); 610 } 611 612 #undef __FUNCT__ 613 #define __FUNCT__ "PCMGGetLevels" 614 /*@ 615 PCMGGetLevels - Gets the number of levels to use with MG. 616 617 Not Collective 618 619 Input Parameter: 620 . pc - the preconditioner context 621 622 Output parameter: 623 . levels - the number of levels 624 625 Level: advanced 626 627 .keywords: MG, get, levels, multigrid 628 629 .seealso: PCMGSetLevels() 630 @*/ 631 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels) 632 { 633 PC_MG **mg; 634 635 PetscFunctionBegin; 636 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 637 PetscValidIntPointer(levels,2); 638 639 mg = (PC_MG**)pc->data; 640 *levels = mg[0]->levels; 641 PetscFunctionReturn(0); 642 } 643 644 #undef __FUNCT__ 645 #define __FUNCT__ "PCMGSetType" 646 /*@ 647 PCMGSetType - Determines the form of multigrid to use: 648 multiplicative, additive, full, or the Kaskade algorithm. 649 650 Collective on PC 651 652 Input Parameters: 653 + pc - the preconditioner context 654 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 655 PC_MG_FULL, PC_MG_KASKADE 656 657 Options Database Key: 658 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 659 additive, full, kaskade 660 661 Level: advanced 662 663 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 664 665 .seealso: PCMGSetLevels() 666 @*/ 667 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form) 668 { 669 PC_MG **mg; 670 671 PetscFunctionBegin; 672 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 673 mg = (PC_MG**)pc->data; 674 675 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 676 mg[0]->am = form; 677 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 678 else pc->ops->applyrichardson = 0; 679 PetscFunctionReturn(0); 680 } 681 682 #undef __FUNCT__ 683 #define __FUNCT__ "PCMGSetCycleType" 684 /*@ 685 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 686 complicated cycling. 687 688 Collective on PC 689 690 Input Parameters: 691 + pc - the multigrid context 692 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 693 694 Options Database Key: 695 $ -pc_mg_cycle_type v or w 696 697 Level: advanced 698 699 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 700 701 .seealso: PCMGSetCycleTypeOnLevel() 702 @*/ 703 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n) 704 { 705 PC_MG **mg; 706 PetscInt i,levels; 707 708 PetscFunctionBegin; 709 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 710 mg = (PC_MG**)pc->data; 711 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 712 levels = mg[0]->levels; 713 714 for (i=0; i<levels; i++) { 715 mg[i]->cycles = n; 716 } 717 PetscFunctionReturn(0); 718 } 719 720 #undef __FUNCT__ 721 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 722 /*@ 723 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 724 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 725 726 Collective on PC 727 728 Input Parameters: 729 + pc - the multigrid context 730 - n - number of cycles (default is 1) 731 732 Options Database Key: 733 $ -pc_mg_multiplicative_cycles n 734 735 Level: advanced 736 737 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 738 739 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 740 741 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 742 @*/ 743 PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 744 { 745 PC_MG **mg; 746 PetscInt i,levels; 747 748 PetscFunctionBegin; 749 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 750 mg = (PC_MG**)pc->data; 751 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 752 levels = mg[0]->levels; 753 754 for (i=0; i<levels; i++) { 755 mg[i]->cyclesperpcapply = n; 756 } 757 PetscFunctionReturn(0); 758 } 759 760 #undef __FUNCT__ 761 #define __FUNCT__ "PCMGSetGalerkin" 762 /*@ 763 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 764 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 765 766 Collective on PC 767 768 Input Parameters: 769 . pc - the multigrid context 770 771 Options Database Key: 772 $ -pc_mg_galerkin 773 774 Level: intermediate 775 776 .keywords: MG, set, Galerkin 777 778 .seealso: PCMGGetGalerkin() 779 780 @*/ 781 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 782 { 783 PC_MG **mg; 784 PetscInt i,levels; 785 786 PetscFunctionBegin; 787 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 788 mg = (PC_MG**)pc->data; 789 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 790 levels = mg[0]->levels; 791 792 for (i=0; i<levels; i++) { 793 mg[i]->galerkin = PETSC_TRUE; 794 } 795 PetscFunctionReturn(0); 796 } 797 798 #undef __FUNCT__ 799 #define __FUNCT__ "PCMGGetGalerkin" 800 /*@ 801 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 802 A_i-1 = r_i * A_i * r_i^t 803 804 Not Collective 805 806 Input Parameter: 807 . pc - the multigrid context 808 809 Output Parameter: 810 . gelerkin - PETSC_TRUE or PETSC_FALSE 811 812 Options Database Key: 813 $ -pc_mg_galerkin 814 815 Level: intermediate 816 817 .keywords: MG, set, Galerkin 818 819 .seealso: PCMGSetGalerkin() 820 821 @*/ 822 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 823 { 824 PC_MG **mg; 825 826 PetscFunctionBegin; 827 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 828 mg = (PC_MG**)pc->data; 829 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 830 *galerkin = mg[0]->galerkin; 831 PetscFunctionReturn(0); 832 } 833 834 #undef __FUNCT__ 835 #define __FUNCT__ "PCMGSetNumberSmoothDown" 836 /*@ 837 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 838 use on all levels. Use PCMGGetSmootherDown() to set different 839 pre-smoothing steps on different levels. 840 841 Collective on PC 842 843 Input Parameters: 844 + mg - the multigrid context 845 - n - the number of smoothing steps 846 847 Options Database Key: 848 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 849 850 Level: advanced 851 852 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 853 854 .seealso: PCMGSetNumberSmoothUp() 855 @*/ 856 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n) 857 { 858 PC_MG **mg; 859 PetscErrorCode ierr; 860 PetscInt i,levels; 861 862 PetscFunctionBegin; 863 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 864 mg = (PC_MG**)pc->data; 865 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 866 levels = mg[0]->levels; 867 868 for (i=1; i<levels; i++) { 869 /* make sure smoother up and down are different */ 870 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 871 ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 872 mg[i]->default_smoothd = n; 873 } 874 PetscFunctionReturn(0); 875 } 876 877 #undef __FUNCT__ 878 #define __FUNCT__ "PCMGSetNumberSmoothUp" 879 /*@ 880 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 881 on all levels. Use PCMGGetSmootherUp() to set different numbers of 882 post-smoothing steps on different levels. 883 884 Collective on PC 885 886 Input Parameters: 887 + mg - the multigrid context 888 - n - the number of smoothing steps 889 890 Options Database Key: 891 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 892 893 Level: advanced 894 895 Note: this does not set a value on the coarsest grid, since we assume that 896 there is no separate smooth up on the coarsest grid. 897 898 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 899 900 .seealso: PCMGSetNumberSmoothDown() 901 @*/ 902 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 903 { 904 PC_MG **mg; 905 PetscErrorCode ierr; 906 PetscInt i,levels; 907 908 PetscFunctionBegin; 909 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 910 mg = (PC_MG**)pc->data; 911 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 912 levels = mg[0]->levels; 913 914 for (i=1; i<levels; i++) { 915 /* make sure smoother up and down are different */ 916 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 917 ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 918 mg[i]->default_smoothu = n; 919 } 920 PetscFunctionReturn(0); 921 } 922 923 /* ----------------------------------------------------------------------------------------*/ 924 925 /*MC 926 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 927 information about the coarser grid matrices and restriction/interpolation operators. 928 929 Options Database Keys: 930 + -pc_mg_levels <nlevels> - number of levels including finest 931 . -pc_mg_cycles v or w 932 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 933 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 934 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 935 . -pc_mg_log - log information about time spent on each level of the solver 936 . -pc_mg_monitor - print information on the multigrid convergence 937 . -pc_mg_galerkin - use Galerkin process to compute coarser operators 938 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 939 to the Socket viewer for reading from Matlab. 940 941 Notes: 942 943 Level: intermediate 944 945 Concepts: multigrid/multilevel 946 947 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 948 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 949 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 950 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 951 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 952 M*/ 953 954 EXTERN_C_BEGIN 955 #undef __FUNCT__ 956 #define __FUNCT__ "PCCreate_MG" 957 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 958 { 959 PetscFunctionBegin; 960 pc->ops->apply = PCApply_MG; 961 pc->ops->setup = PCSetUp_MG; 962 pc->ops->destroy = PCDestroy_MG; 963 pc->ops->setfromoptions = PCSetFromOptions_MG; 964 pc->ops->view = PCView_MG; 965 966 pc->data = (void*)0; 967 PetscFunctionReturn(0); 968 } 969 EXTERN_C_END 970