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