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