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