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