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