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