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