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