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 #if defined(PETSC_USE_LOG) 418 { 419 const char *sname = "MG Apply"; 420 PetscStageLog stageLog; 421 PetscInt st; 422 423 PetscFunctionBegin; 424 ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 425 for (st = 0; st < stageLog->numStages; ++st) { 426 PetscBool same; 427 428 ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 429 if (same) {mg->stageApply = st;} 430 } 431 if (!mg->stageApply) { 432 ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 433 } 434 } 435 #endif 436 } 437 ierr = PetscOptionsTail();CHKERRQ(ierr); 438 PetscFunctionReturn(0); 439 } 440 441 const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 442 const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 443 444 #undef __FUNCT__ 445 #define __FUNCT__ "PCView_MG" 446 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 447 { 448 PC_MG *mg = (PC_MG*)pc->data; 449 PC_MG_Levels **mglevels = mg->levels; 450 PetscErrorCode ierr; 451 PetscInt levels = mglevels[0]->levels,i; 452 PetscBool iascii,isbinary,isdraw; 453 454 PetscFunctionBegin; 455 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 456 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 457 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 458 if (iascii) { 459 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); 460 if (mg->am == PC_MG_MULTIPLICATIVE) { 461 ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 462 } 463 if (mg->galerkin) { 464 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 465 } else { 466 ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 467 } 468 for (i=0; i<levels; i++) { 469 if (!i) { 470 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 471 } else { 472 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 473 } 474 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 475 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 476 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 477 if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 478 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 479 } else if (i){ 480 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 481 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 482 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 483 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 484 } 485 } 486 } else if (isbinary) { 487 for (i=levels-1; i>=0; i--) { 488 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 489 if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 490 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 491 } 492 } 493 } else if (isdraw) { 494 PetscDraw draw; 495 PetscReal x,w,y,bottom,th; 496 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 497 ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 498 ierr = PetscDrawStringGetSize(draw,PETSC_NULL,&th);CHKERRQ(ierr); 499 bottom = y - th; 500 for (i=levels-1; i>=0; i--) { 501 if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 502 ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 503 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 504 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 505 } else { 506 w = 0.5*PetscMin(1.0-x,x); 507 ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 508 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 509 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 510 ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 511 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 512 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 513 } 514 ierr = PetscDrawGetBoundingBox(draw,PETSC_NULL,&bottom,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 515 bottom -= th; 516 } 517 } 518 PetscFunctionReturn(0); 519 } 520 521 #include <petsc-private/dmimpl.h> 522 #include <petsc-private/kspimpl.h> 523 524 /* 525 Calls setup for the KSP on each level 526 */ 527 #undef __FUNCT__ 528 #define __FUNCT__ "PCSetUp_MG" 529 PetscErrorCode PCSetUp_MG(PC pc) 530 { 531 PC_MG *mg = (PC_MG*)pc->data; 532 PC_MG_Levels **mglevels = mg->levels; 533 PetscErrorCode ierr; 534 PetscInt i,n = mglevels[0]->levels; 535 PC cpc; 536 PetscBool preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset; 537 Mat dA,dB; 538 MatStructure uflag; 539 Vec tvec; 540 DM *dms; 541 PetscViewer viewer = 0; 542 543 PetscFunctionBegin; 544 /* FIX: Move this to PCSetFromOptions_MG? */ 545 if (mg->usedmfornumberoflevels) { 546 PetscInt levels; 547 ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 548 levels++; 549 if (levels > n) { /* the problem is now being solved on a finer grid */ 550 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 551 n = levels; 552 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 553 mglevels = mg->levels; 554 } 555 } 556 ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 557 558 559 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 560 /* so use those from global PC */ 561 /* Is this what we always want? What if user wants to keep old one? */ 562 ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 563 if (opsset) { 564 Mat mmat; 565 ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr); 566 if (mmat == pc->pmat) opsset = PETSC_FALSE; 567 } 568 if (!opsset) { 569 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); 570 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 571 } 572 573 /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */ 574 if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) { 575 /* construct the interpolation from the DMs */ 576 Mat p; 577 Vec rscale; 578 ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 579 dms[n-1] = pc->dm; 580 for (i=n-2; i>-1; i--) { 581 DMKSP kdm; 582 ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr); 583 ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 584 if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 585 ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 586 /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 587 * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 588 kdm->ops->computerhs = PETSC_NULL; 589 kdm->rhsctx = PETSC_NULL; 590 if (!mglevels[i+1]->interpolate) { 591 ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 592 ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 593 if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 594 ierr = VecDestroy(&rscale);CHKERRQ(ierr); 595 ierr = MatDestroy(&p);CHKERRQ(ierr); 596 } 597 } 598 599 for (i=n-2; i>-1; i--) { 600 ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 601 } 602 ierr = PetscFree(dms);CHKERRQ(ierr); 603 } 604 605 if (pc->dm && !pc->setupcalled) { 606 /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */ 607 ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 608 ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 609 } 610 611 if (mg->galerkin == 1) { 612 Mat B; 613 /* currently only handle case where mat and pmat are the same on coarser levels */ 614 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 615 if (!pc->setupcalled) { 616 for (i=n-2; i>-1; i--) { 617 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 618 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 619 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 620 dB = B; 621 } 622 if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 623 } else { 624 for (i=n-2; i>-1; i--) { 625 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 626 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 627 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 628 dB = B; 629 } 630 } 631 } else if (!mg->galerkin && pc->dm && pc->dm->x) { 632 /* need to restrict Jacobian location to coarser meshes for evaluation */ 633 for (i=n-2;i>-1; i--) { 634 Mat R; 635 Vec rscale; 636 if (!mglevels[i]->smoothd->dm->x) { 637 Vec *vecs; 638 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr); 639 mglevels[i]->smoothd->dm->x = vecs[0]; 640 ierr = PetscFree(vecs);CHKERRQ(ierr); 641 } 642 ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 643 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 644 ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 645 ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 646 } 647 } 648 if (!mg->galerkin && pc->dm) { 649 for (i=n-2;i>=0; i--) { 650 DM dmfine,dmcoarse; 651 Mat Restrict,Inject; 652 Vec rscale; 653 ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 654 ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 655 ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 656 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 657 Inject = PETSC_NULL; /* Callback should create it if it needs Injection */ 658 ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 659 } 660 } 661 662 if (!pc->setupcalled) { 663 for (i=0; i<n; i++) { 664 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 665 } 666 for (i=1; i<n; i++) { 667 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 668 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 669 } 670 } 671 for (i=1; i<n; i++) { 672 ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 673 ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 674 } 675 for (i=0; i<n-1; i++) { 676 if (!mglevels[i]->b) { 677 Vec *vec; 678 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 679 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 680 ierr = VecDestroy(vec);CHKERRQ(ierr); 681 ierr = PetscFree(vec);CHKERRQ(ierr); 682 } 683 if (!mglevels[i]->r && i) { 684 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 685 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 686 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 687 } 688 if (!mglevels[i]->x) { 689 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 690 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 691 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 692 } 693 } 694 if (n != 1 && !mglevels[n-1]->r) { 695 /* PCMGSetR() on the finest level if user did not supply it */ 696 Vec *vec; 697 ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 698 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 699 ierr = VecDestroy(vec);CHKERRQ(ierr); 700 ierr = PetscFree(vec);CHKERRQ(ierr); 701 } 702 } 703 704 if (pc->dm) { 705 /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 706 for (i=0; i<n-1; i++){ 707 if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 708 } 709 } 710 711 for (i=1; i<n; i++) { 712 if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 713 /* if doing only down then initial guess is zero */ 714 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 715 } 716 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 717 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 718 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 719 if (!mglevels[i]->residual) { 720 Mat mat; 721 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 722 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 723 } 724 } 725 for (i=1; i<n; i++) { 726 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 727 Mat downmat,downpmat; 728 MatStructure matflag; 729 730 /* check if operators have been set for up, if not use down operators to set them */ 731 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 732 if (!opsset) { 733 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 734 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 735 } 736 737 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 738 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 739 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 740 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 741 } 742 } 743 744 /* 745 If coarse solver is not direct method then DO NOT USE preonly 746 */ 747 ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 748 if (preonly) { 749 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 750 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 751 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 752 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 753 if (!lu && !redundant && !cholesky && !svd) { 754 ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 755 } 756 } 757 758 if (!pc->setupcalled) { 759 ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 760 } 761 762 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 763 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 764 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 765 766 /* 767 Dump the interpolation/restriction matrices plus the 768 Jacobian/stiffness on each level. This allows MATLAB users to 769 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 770 771 Only support one or the other at the same time. 772 */ 773 #if defined(PETSC_USE_SOCKET_VIEWER) 774 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 775 if (dump) { 776 viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 777 } 778 dump = PETSC_FALSE; 779 #endif 780 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 781 if (dump) { 782 783 viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 784 } 785 786 if (viewer) { 787 for (i=1; i<n; i++) { 788 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 789 } 790 for (i=0; i<n; i++) { 791 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 792 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 793 } 794 } 795 PetscFunctionReturn(0); 796 } 797 798 /* -------------------------------------------------------------------------------------*/ 799 800 #undef __FUNCT__ 801 #define __FUNCT__ "PCMGGetLevels" 802 /*@ 803 PCMGGetLevels - Gets the number of levels to use with MG. 804 805 Not Collective 806 807 Input Parameter: 808 . pc - the preconditioner context 809 810 Output parameter: 811 . levels - the number of levels 812 813 Level: advanced 814 815 .keywords: MG, get, levels, multigrid 816 817 .seealso: PCMGSetLevels() 818 @*/ 819 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 820 { 821 PC_MG *mg = (PC_MG*)pc->data; 822 823 PetscFunctionBegin; 824 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 825 PetscValidIntPointer(levels,2); 826 *levels = mg->nlevels; 827 PetscFunctionReturn(0); 828 } 829 830 #undef __FUNCT__ 831 #define __FUNCT__ "PCMGSetType" 832 /*@ 833 PCMGSetType - Determines the form of multigrid to use: 834 multiplicative, additive, full, or the Kaskade algorithm. 835 836 Logically Collective on PC 837 838 Input Parameters: 839 + pc - the preconditioner context 840 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 841 PC_MG_FULL, PC_MG_KASKADE 842 843 Options Database Key: 844 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 845 additive, full, kaskade 846 847 Level: advanced 848 849 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 850 851 .seealso: PCMGSetLevels() 852 @*/ 853 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 854 { 855 PC_MG *mg = (PC_MG*)pc->data; 856 857 PetscFunctionBegin; 858 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 859 PetscValidLogicalCollectiveEnum(pc,form,2); 860 mg->am = form; 861 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 862 else pc->ops->applyrichardson = 0; 863 PetscFunctionReturn(0); 864 } 865 866 #undef __FUNCT__ 867 #define __FUNCT__ "PCMGSetCycleType" 868 /*@ 869 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 870 complicated cycling. 871 872 Logically Collective on PC 873 874 Input Parameters: 875 + pc - the multigrid context 876 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 877 878 Options Database Key: 879 $ -pc_mg_cycle_type v or w 880 881 Level: advanced 882 883 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 884 885 .seealso: PCMGSetCycleTypeOnLevel() 886 @*/ 887 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 888 { 889 PC_MG *mg = (PC_MG*)pc->data; 890 PC_MG_Levels **mglevels = mg->levels; 891 PetscInt i,levels; 892 893 PetscFunctionBegin; 894 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 895 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 896 PetscValidLogicalCollectiveInt(pc,n,2); 897 levels = mglevels[0]->levels; 898 899 for (i=0; i<levels; i++) { 900 mglevels[i]->cycles = n; 901 } 902 PetscFunctionReturn(0); 903 } 904 905 #undef __FUNCT__ 906 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 907 /*@ 908 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 909 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 910 911 Logically Collective on PC 912 913 Input Parameters: 914 + pc - the multigrid context 915 - n - number of cycles (default is 1) 916 917 Options Database Key: 918 $ -pc_mg_multiplicative_cycles n 919 920 Level: advanced 921 922 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 923 924 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 925 926 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 927 @*/ 928 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 929 { 930 PC_MG *mg = (PC_MG*)pc->data; 931 PC_MG_Levels **mglevels = mg->levels; 932 PetscInt i,levels; 933 934 PetscFunctionBegin; 935 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 936 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 937 PetscValidLogicalCollectiveInt(pc,n,2); 938 levels = mglevels[0]->levels; 939 940 for (i=0; i<levels; i++) { 941 mg->cyclesperpcapply = n; 942 } 943 PetscFunctionReturn(0); 944 } 945 946 #undef __FUNCT__ 947 #define __FUNCT__ "PCMGSetGalerkin" 948 /*@ 949 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 950 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 951 952 Logically Collective on PC 953 954 Input Parameters: 955 + pc - the multigrid context 956 - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 957 958 Options Database Key: 959 $ -pc_mg_galerkin 960 961 Level: intermediate 962 963 .keywords: MG, set, Galerkin 964 965 .seealso: PCMGGetGalerkin() 966 967 @*/ 968 PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 969 { 970 PC_MG *mg = (PC_MG*)pc->data; 971 972 PetscFunctionBegin; 973 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 974 mg->galerkin = use ? 1 : 0; 975 PetscFunctionReturn(0); 976 } 977 978 #undef __FUNCT__ 979 #define __FUNCT__ "PCMGGetGalerkin" 980 /*@ 981 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 982 A_i-1 = r_i * A_i * r_i^t 983 984 Not Collective 985 986 Input Parameter: 987 . pc - the multigrid context 988 989 Output Parameter: 990 . gelerkin - PETSC_TRUE or PETSC_FALSE 991 992 Options Database Key: 993 $ -pc_mg_galerkin 994 995 Level: intermediate 996 997 .keywords: MG, set, Galerkin 998 999 .seealso: PCMGSetGalerkin() 1000 1001 @*/ 1002 PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 1003 { 1004 PC_MG *mg = (PC_MG*)pc->data; 1005 1006 PetscFunctionBegin; 1007 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1008 *galerkin = (PetscBool)mg->galerkin; 1009 PetscFunctionReturn(0); 1010 } 1011 1012 #undef __FUNCT__ 1013 #define __FUNCT__ "PCMGSetNumberSmoothDown" 1014 /*@ 1015 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 1016 use on all levels. Use PCMGGetSmootherDown() to set different 1017 pre-smoothing steps on different levels. 1018 1019 Logically Collective on PC 1020 1021 Input Parameters: 1022 + mg - the multigrid context 1023 - n - the number of smoothing steps 1024 1025 Options Database Key: 1026 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 1027 1028 Level: advanced 1029 1030 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 1031 1032 .seealso: PCMGSetNumberSmoothUp() 1033 @*/ 1034 PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 1035 { 1036 PC_MG *mg = (PC_MG*)pc->data; 1037 PC_MG_Levels **mglevels = mg->levels; 1038 PetscErrorCode ierr; 1039 PetscInt i,levels; 1040 1041 PetscFunctionBegin; 1042 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1043 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1044 PetscValidLogicalCollectiveInt(pc,n,2); 1045 levels = mglevels[0]->levels; 1046 1047 for (i=1; i<levels; i++) { 1048 /* make sure smoother up and down are different */ 1049 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1050 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1051 mg->default_smoothd = n; 1052 } 1053 PetscFunctionReturn(0); 1054 } 1055 1056 #undef __FUNCT__ 1057 #define __FUNCT__ "PCMGSetNumberSmoothUp" 1058 /*@ 1059 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 1060 on all levels. Use PCMGGetSmootherUp() to set different numbers of 1061 post-smoothing steps on different levels. 1062 1063 Logically Collective on PC 1064 1065 Input Parameters: 1066 + mg - the multigrid context 1067 - n - the number of smoothing steps 1068 1069 Options Database Key: 1070 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 1071 1072 Level: advanced 1073 1074 Note: this does not set a value on the coarsest grid, since we assume that 1075 there is no separate smooth up on the coarsest grid. 1076 1077 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1078 1079 .seealso: PCMGSetNumberSmoothDown() 1080 @*/ 1081 PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 1082 { 1083 PC_MG *mg = (PC_MG*)pc->data; 1084 PC_MG_Levels **mglevels = mg->levels; 1085 PetscErrorCode ierr; 1086 PetscInt i,levels; 1087 1088 PetscFunctionBegin; 1089 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1090 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1091 PetscValidLogicalCollectiveInt(pc,n,2); 1092 levels = mglevels[0]->levels; 1093 1094 for (i=1; i<levels; i++) { 1095 /* make sure smoother up and down are different */ 1096 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1097 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1098 mg->default_smoothu = n; 1099 } 1100 PetscFunctionReturn(0); 1101 } 1102 1103 /* ----------------------------------------------------------------------------------------*/ 1104 1105 /*MC 1106 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1107 information about the coarser grid matrices and restriction/interpolation operators. 1108 1109 Options Database Keys: 1110 + -pc_mg_levels <nlevels> - number of levels including finest 1111 . -pc_mg_cycles v or w 1112 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 1113 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 1114 . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 1115 . -pc_mg_log - log information about time spent on each level of the solver 1116 . -pc_mg_monitor - print information on the multigrid convergence 1117 . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 1118 . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 1119 . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1120 to the Socket viewer for reading from MATLAB. 1121 - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 1122 to the binary output file called binaryoutput 1123 1124 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 1125 1126 When run with a single level the smoother options are used on that level NOT the coarse grid solver options 1127 1128 Level: intermediate 1129 1130 Concepts: multigrid/multilevel 1131 1132 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1133 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1134 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1135 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1136 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1137 M*/ 1138 1139 EXTERN_C_BEGIN 1140 #undef __FUNCT__ 1141 #define __FUNCT__ "PCCreate_MG" 1142 PetscErrorCode PCCreate_MG(PC pc) 1143 { 1144 PC_MG *mg; 1145 PetscErrorCode ierr; 1146 1147 PetscFunctionBegin; 1148 ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1149 pc->data = (void*)mg; 1150 mg->nlevels = -1; 1151 1152 pc->ops->apply = PCApply_MG; 1153 pc->ops->setup = PCSetUp_MG; 1154 pc->ops->reset = PCReset_MG; 1155 pc->ops->destroy = PCDestroy_MG; 1156 pc->ops->setfromoptions = PCSetFromOptions_MG; 1157 pc->ops->view = PCView_MG; 1158 PetscFunctionReturn(0); 1159 } 1160 EXTERN_C_END 1161