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 if (pc->dm && !pc->setupcalled) { 509 /* construct the interpolation from the DMs */ 510 Mat p; 511 Vec rscale; 512 ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 513 dms[n-1] = pc->dm; 514 for (i=n-2; i>-1; i--) { 515 ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr); 516 ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 517 if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 518 ierr = DMSetFunction(dms[i],0); 519 ierr = DMSetInitialGuess(dms[i],0); 520 if (!mglevels[i+1]->interpolate) { 521 ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 522 ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 523 if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 524 ierr = VecDestroy(&rscale);CHKERRQ(ierr); 525 ierr = MatDestroy(&p);CHKERRQ(ierr); 526 } 527 } 528 529 for (i=n-2; i>-1; i--) { 530 ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 531 } 532 ierr = PetscFree(dms);CHKERRQ(ierr); 533 534 /* finest smoother also gets DM but it is not active */ 535 ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 536 ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 537 } 538 539 if (mg->galerkin == 1) { 540 Mat B; 541 /* currently only handle case where mat and pmat are the same on coarser levels */ 542 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 543 if (!pc->setupcalled) { 544 for (i=n-2; i>-1; i--) { 545 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 546 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 547 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 548 dB = B; 549 } 550 if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 551 } else { 552 for (i=n-2; i>-1; i--) { 553 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 554 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 555 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 556 dB = B; 557 } 558 } 559 } else if (pc->dm && pc->dm->x) { 560 /* need to restrict Jacobian location to coarser meshes for evaluation */ 561 for (i=n-2;i>-1; i--) { 562 Mat R; 563 Vec rscale; 564 if (!mglevels[i]->smoothd->dm->x) { 565 Vec *vecs; 566 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr); 567 mglevels[i]->smoothd->dm->x = vecs[0]; 568 ierr = PetscFree(vecs);CHKERRQ(ierr); 569 } 570 ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 571 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 572 ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 573 ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 574 } 575 } 576 577 if (!pc->setupcalled) { 578 for (i=0; i<n; i++) { 579 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 580 } 581 for (i=1; i<n; i++) { 582 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 583 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 584 } 585 } 586 for (i=1; i<n; i++) { 587 ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 588 ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 589 } 590 for (i=0; i<n-1; i++) { 591 if (!mglevels[i]->b) { 592 Vec *vec; 593 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 594 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 595 ierr = VecDestroy(vec);CHKERRQ(ierr); 596 ierr = PetscFree(vec);CHKERRQ(ierr); 597 } 598 if (!mglevels[i]->r && i) { 599 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 600 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 601 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 602 } 603 if (!mglevels[i]->x) { 604 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 605 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 606 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 607 } 608 } 609 if (n != 1 && !mglevels[n-1]->r) { 610 /* PCMGSetR() on the finest level if user did not supply it */ 611 Vec *vec; 612 ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 613 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 614 ierr = VecDestroy(vec);CHKERRQ(ierr); 615 ierr = PetscFree(vec);CHKERRQ(ierr); 616 } 617 } 618 619 620 for (i=1; i<n; i++) { 621 if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 622 /* if doing only down then initial guess is zero */ 623 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 624 } 625 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 626 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 627 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 628 if (!mglevels[i]->residual) { 629 Mat mat; 630 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 631 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 632 } 633 } 634 for (i=1; i<n; i++) { 635 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 636 Mat downmat,downpmat; 637 MatStructure matflag; 638 PetscBool opsset; 639 640 /* check if operators have been set for up, if not use down operators to set them */ 641 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 642 if (!opsset) { 643 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 644 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 645 } 646 647 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 648 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 649 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 650 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 651 } 652 } 653 654 /* 655 If coarse solver is not direct method then DO NOT USE preonly 656 */ 657 ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 658 if (preonly) { 659 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 660 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 661 ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 662 ierr = PetscTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 663 if (!lu && !redundant && !cholesky && !svd) { 664 ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 665 } 666 } 667 668 if (!pc->setupcalled) { 669 ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 670 } 671 672 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 673 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 674 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 675 676 /* 677 Dump the interpolation/restriction matrices plus the 678 Jacobian/stiffness on each level. This allows MATLAB users to 679 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 680 681 Only support one or the other at the same time. 682 */ 683 #if defined(PETSC_USE_SOCKET_VIEWER) 684 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 685 if (dump) { 686 viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 687 } 688 dump = PETSC_FALSE; 689 #endif 690 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 691 if (dump) { 692 viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 693 } 694 695 if (viewer) { 696 for (i=1; i<n; i++) { 697 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 698 } 699 for (i=0; i<n; i++) { 700 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 701 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 702 } 703 } 704 PetscFunctionReturn(0); 705 } 706 707 /* -------------------------------------------------------------------------------------*/ 708 709 #undef __FUNCT__ 710 #define __FUNCT__ "PCMGGetLevels" 711 /*@ 712 PCMGGetLevels - Gets the number of levels to use with MG. 713 714 Not Collective 715 716 Input Parameter: 717 . pc - the preconditioner context 718 719 Output parameter: 720 . levels - the number of levels 721 722 Level: advanced 723 724 .keywords: MG, get, levels, multigrid 725 726 .seealso: PCMGSetLevels() 727 @*/ 728 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 729 { 730 PC_MG *mg = (PC_MG*)pc->data; 731 732 PetscFunctionBegin; 733 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 734 PetscValidIntPointer(levels,2); 735 *levels = mg->nlevels; 736 PetscFunctionReturn(0); 737 } 738 739 #undef __FUNCT__ 740 #define __FUNCT__ "PCMGSetType" 741 /*@ 742 PCMGSetType - Determines the form of multigrid to use: 743 multiplicative, additive, full, or the Kaskade algorithm. 744 745 Logically Collective on PC 746 747 Input Parameters: 748 + pc - the preconditioner context 749 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 750 PC_MG_FULL, PC_MG_KASKADE 751 752 Options Database Key: 753 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 754 additive, full, kaskade 755 756 Level: advanced 757 758 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 759 760 .seealso: PCMGSetLevels() 761 @*/ 762 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 763 { 764 PC_MG *mg = (PC_MG*)pc->data; 765 766 PetscFunctionBegin; 767 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 768 PetscValidLogicalCollectiveEnum(pc,form,2); 769 mg->am = form; 770 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 771 else pc->ops->applyrichardson = 0; 772 PetscFunctionReturn(0); 773 } 774 775 #undef __FUNCT__ 776 #define __FUNCT__ "PCMGSetCycleType" 777 /*@ 778 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 779 complicated cycling. 780 781 Logically Collective on PC 782 783 Input Parameters: 784 + pc - the multigrid context 785 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 786 787 Options Database Key: 788 $ -pc_mg_cycle_type v or w 789 790 Level: advanced 791 792 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 793 794 .seealso: PCMGSetCycleTypeOnLevel() 795 @*/ 796 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 797 { 798 PC_MG *mg = (PC_MG*)pc->data; 799 PC_MG_Levels **mglevels = mg->levels; 800 PetscInt i,levels; 801 802 PetscFunctionBegin; 803 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 804 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 805 PetscValidLogicalCollectiveInt(pc,n,2); 806 levels = mglevels[0]->levels; 807 808 for (i=0; i<levels; i++) { 809 mglevels[i]->cycles = n; 810 } 811 PetscFunctionReturn(0); 812 } 813 814 #undef __FUNCT__ 815 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 816 /*@ 817 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 818 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 819 820 Logically Collective on PC 821 822 Input Parameters: 823 + pc - the multigrid context 824 - n - number of cycles (default is 1) 825 826 Options Database Key: 827 $ -pc_mg_multiplicative_cycles n 828 829 Level: advanced 830 831 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 832 833 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 834 835 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 836 @*/ 837 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 838 { 839 PC_MG *mg = (PC_MG*)pc->data; 840 PC_MG_Levels **mglevels = mg->levels; 841 PetscInt i,levels; 842 843 PetscFunctionBegin; 844 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 845 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 846 PetscValidLogicalCollectiveInt(pc,n,2); 847 levels = mglevels[0]->levels; 848 849 for (i=0; i<levels; i++) { 850 mg->cyclesperpcapply = n; 851 } 852 PetscFunctionReturn(0); 853 } 854 855 #undef __FUNCT__ 856 #define __FUNCT__ "PCMGSetGalerkin" 857 /*@ 858 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 859 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 860 861 Logically Collective on PC 862 863 Input Parameters: 864 + pc - the multigrid context 865 - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 866 867 Options Database Key: 868 $ -pc_mg_galerkin 869 870 Level: intermediate 871 872 .keywords: MG, set, Galerkin 873 874 .seealso: PCMGGetGalerkin() 875 876 @*/ 877 PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 878 { 879 PC_MG *mg = (PC_MG*)pc->data; 880 881 PetscFunctionBegin; 882 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 883 mg->galerkin = (PetscInt)use; 884 PetscFunctionReturn(0); 885 } 886 887 #undef __FUNCT__ 888 #define __FUNCT__ "PCMGGetGalerkin" 889 /*@ 890 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 891 A_i-1 = r_i * A_i * r_i^t 892 893 Not Collective 894 895 Input Parameter: 896 . pc - the multigrid context 897 898 Output Parameter: 899 . gelerkin - PETSC_TRUE or PETSC_FALSE 900 901 Options Database Key: 902 $ -pc_mg_galerkin 903 904 Level: intermediate 905 906 .keywords: MG, set, Galerkin 907 908 .seealso: PCMGSetGalerkin() 909 910 @*/ 911 PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 912 { 913 PC_MG *mg = (PC_MG*)pc->data; 914 915 PetscFunctionBegin; 916 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 917 *galerkin = (PetscBool)mg->galerkin; 918 PetscFunctionReturn(0); 919 } 920 921 #undef __FUNCT__ 922 #define __FUNCT__ "PCMGSetNumberSmoothDown" 923 /*@ 924 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 925 use on all levels. Use PCMGGetSmootherDown() to set different 926 pre-smoothing steps on different levels. 927 928 Logically Collective on PC 929 930 Input Parameters: 931 + mg - the multigrid context 932 - n - the number of smoothing steps 933 934 Options Database Key: 935 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 936 937 Level: advanced 938 939 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 940 941 .seealso: PCMGSetNumberSmoothUp() 942 @*/ 943 PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 944 { 945 PC_MG *mg = (PC_MG*)pc->data; 946 PC_MG_Levels **mglevels = mg->levels; 947 PetscErrorCode ierr; 948 PetscInt i,levels; 949 950 PetscFunctionBegin; 951 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 952 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 953 PetscValidLogicalCollectiveInt(pc,n,2); 954 levels = mglevels[0]->levels; 955 956 for (i=1; i<levels; i++) { 957 /* make sure smoother up and down are different */ 958 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 959 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 960 mg->default_smoothd = n; 961 } 962 PetscFunctionReturn(0); 963 } 964 965 #undef __FUNCT__ 966 #define __FUNCT__ "PCMGSetNumberSmoothUp" 967 /*@ 968 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 969 on all levels. Use PCMGGetSmootherUp() to set different numbers of 970 post-smoothing steps on different levels. 971 972 Logically Collective on PC 973 974 Input Parameters: 975 + mg - the multigrid context 976 - n - the number of smoothing steps 977 978 Options Database Key: 979 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 980 981 Level: advanced 982 983 Note: this does not set a value on the coarsest grid, since we assume that 984 there is no separate smooth up on the coarsest grid. 985 986 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 987 988 .seealso: PCMGSetNumberSmoothDown() 989 @*/ 990 PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 991 { 992 PC_MG *mg = (PC_MG*)pc->data; 993 PC_MG_Levels **mglevels = mg->levels; 994 PetscErrorCode ierr; 995 PetscInt i,levels; 996 997 PetscFunctionBegin; 998 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 999 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1000 PetscValidLogicalCollectiveInt(pc,n,2); 1001 levels = mglevels[0]->levels; 1002 1003 for (i=1; i<levels; i++) { 1004 /* make sure smoother up and down are different */ 1005 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1006 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1007 mg->default_smoothu = n; 1008 } 1009 PetscFunctionReturn(0); 1010 } 1011 1012 /* ----------------------------------------------------------------------------------------*/ 1013 1014 /*MC 1015 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1016 information about the coarser grid matrices and restriction/interpolation operators. 1017 1018 Options Database Keys: 1019 + -pc_mg_levels <nlevels> - number of levels including finest 1020 . -pc_mg_cycles v or w 1021 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 1022 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 1023 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 1024 . -pc_mg_log - log information about time spent on each level of the solver 1025 . -pc_mg_monitor - print information on the multigrid convergence 1026 . -pc_mg_galerkin - use Galerkin process to compute coarser operators 1027 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1028 to the Socket viewer for reading from MATLAB. 1029 1030 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 1031 1032 Level: intermediate 1033 1034 Concepts: multigrid/multilevel 1035 1036 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC 1037 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1038 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1039 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1040 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1041 M*/ 1042 1043 EXTERN_C_BEGIN 1044 #undef __FUNCT__ 1045 #define __FUNCT__ "PCCreate_MG" 1046 PetscErrorCode PCCreate_MG(PC pc) 1047 { 1048 PC_MG *mg; 1049 PetscErrorCode ierr; 1050 1051 PetscFunctionBegin; 1052 ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1053 pc->data = (void*)mg; 1054 mg->nlevels = -1; 1055 1056 pc->ops->apply = PCApply_MG; 1057 pc->ops->setup = PCSetUp_MG; 1058 pc->ops->reset = PCReset_MG; 1059 pc->ops->destroy = PCDestroy_MG; 1060 pc->ops->setfromoptions = PCSetFromOptions_MG; 1061 pc->ops->view = PCView_MG; 1062 PetscFunctionReturn(0); 1063 } 1064 EXTERN_C_END 1065