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