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 } 131 132 for (i=0; i<n; i++) { 133 ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 134 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 135 ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 136 } 137 ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 138 } 139 } 140 PetscFunctionReturn(0); 141 } 142 143 #undef __FUNCT__ 144 #define __FUNCT__ "PCMGSetLevels" 145 /*@C 146 PCMGSetLevels - Sets the number of levels to use with MG. 147 Must be called before any other MG routine. 148 149 Logically Collective on PC 150 151 Input Parameters: 152 + pc - the preconditioner context 153 . levels - the number of levels 154 - comms - optional communicators for each level; this is to allow solving the coarser problems 155 on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 156 157 Level: intermediate 158 159 Notes: 160 If the number of levels is one then the multigrid uses the -mg_levels prefix 161 for setting the level options rather than the -mg_coarse prefix. 162 163 .keywords: MG, set, levels, multigrid 164 165 .seealso: PCMGSetType(), PCMGGetLevels() 166 @*/ 167 PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 168 { 169 PetscErrorCode ierr; 170 PC_MG *mg = (PC_MG*)pc->data; 171 MPI_Comm comm = ((PetscObject)pc)->comm; 172 PC_MG_Levels **mglevels = mg->levels; 173 PetscInt i; 174 PetscMPIInt size; 175 const char *prefix; 176 PC ipc; 177 PetscInt n; 178 179 PetscFunctionBegin; 180 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 181 PetscValidLogicalCollectiveInt(pc,levels,2); 182 if (mg->nlevels == levels) PetscFunctionReturn(0); 183 if (mglevels) { 184 /* changing the number of levels so free up the previous stuff */ 185 ierr = PCReset_MG(pc);CHKERRQ(ierr); 186 n = mglevels[0]->levels; 187 for (i=0; i<n; i++) { 188 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 189 ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 190 } 191 ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 192 ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 193 } 194 ierr = PetscFree(mg->levels);CHKERRQ(ierr); 195 } 196 197 mg->nlevels = levels; 198 mg->galerkin = PETSC_FALSE; 199 mg->galerkinused = PETSC_FALSE; 200 201 ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr); 202 ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 203 204 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 205 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 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 PetscFunctionReturn(0); 334 } 335 336 337 #undef __FUNCT__ 338 #define __FUNCT__ "PCSetFromOptions_MG" 339 PetscErrorCode PCSetFromOptions_MG(PC pc) 340 { 341 PetscErrorCode ierr; 342 PetscInt m,levels = 1,cycles; 343 PetscBool flg; 344 PC_MG *mg = (PC_MG*)pc->data; 345 PC_MG_Levels **mglevels = mg->levels; 346 PCMGType mgtype; 347 PCMGCycleType mgctype; 348 349 PetscFunctionBegin; 350 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 351 if (!mg->levels) { 352 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 353 if (!flg && pc->dm) { 354 ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 355 levels++; 356 mg->usedmfornumberoflevels = PETSC_TRUE; 357 } 358 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 359 } 360 mglevels = mg->levels; 361 362 mgctype = (PCMGCycleType) mglevels[0]->cycles; 363 ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 364 if (flg) { 365 ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 366 }; 367 flg = PETSC_FALSE; 368 ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 369 if (flg) { 370 ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 371 } 372 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 373 if (flg) { 374 ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 375 } 376 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 377 if (flg) { 378 ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 379 } 380 mgtype = mg->am; 381 ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 382 if (flg) { 383 ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 384 } 385 if (mg->am == PC_MG_MULTIPLICATIVE) { 386 ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 387 if (flg) { 388 ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 389 } 390 } 391 flg = PETSC_FALSE; 392 ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 393 if (flg) { 394 PetscInt i; 395 char eventname[128]; 396 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 397 levels = mglevels[0]->levels; 398 for (i=0; i<levels; i++) { 399 sprintf(eventname,"MGSetup Level %d",(int)i); 400 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 401 sprintf(eventname,"MGSmooth Level %d",(int)i); 402 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 403 if (i) { 404 sprintf(eventname,"MGResid Level %d",(int)i); 405 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 406 sprintf(eventname,"MGInterp Level %d",(int)i); 407 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 408 } 409 } 410 } 411 ierr = PetscOptionsTail();CHKERRQ(ierr); 412 PetscFunctionReturn(0); 413 } 414 415 const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 416 const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 417 418 #undef __FUNCT__ 419 #define __FUNCT__ "PCView_MG" 420 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 421 { 422 PC_MG *mg = (PC_MG*)pc->data; 423 PC_MG_Levels **mglevels = mg->levels; 424 PetscErrorCode ierr; 425 PetscInt levels = mglevels[0]->levels,i; 426 PetscBool iascii; 427 428 PetscFunctionBegin; 429 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 430 if (iascii) { 431 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); 432 if (mg->am == PC_MG_MULTIPLICATIVE) { 433 ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 434 } 435 if (mg->galerkin) { 436 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 437 } else { 438 ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 439 } 440 for (i=0; i<levels; i++) { 441 if (!i) { 442 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr); 443 } else { 444 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr); 445 } 446 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 447 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 448 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 449 if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 450 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 451 } else if (i){ 452 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr); 453 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 454 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 455 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 456 } 457 } 458 } else { 459 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 460 } 461 PetscFunctionReturn(0); 462 } 463 464 #include <private/dmimpl.h> 465 #include <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,mpc; 479 PetscBool preonly,lu,redundant,cholesky,svd,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset; 480 PetscViewerASCIIMonitor ascii; 481 PetscViewer viewer = PETSC_NULL; 482 MPI_Comm comm; 483 Mat dA,dB; 484 MatStructure uflag; 485 Vec tvec; 486 DM *dms; 487 488 PetscFunctionBegin; 489 if (mg->usedmfornumberoflevels) { 490 PetscInt levels; 491 ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 492 levels++; 493 if (levels > n) { /* the problem is now being solved on a finer grid */ 494 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 495 n = levels; 496 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 497 mglevels = mg->levels; 498 } 499 } 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 ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 507 ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr); 508 if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) { 509 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); 510 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 511 } 512 513 if (pc->dm && !pc->setupcalled) { 514 /* construct the interpolation from the DMs */ 515 Mat p; 516 ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 517 dms[n-1] = pc->dm; 518 for (i=n-2; i>-1; i--) { 519 ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr); 520 ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 521 if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 522 ierr = DMSetFunction(dms[i],0); 523 ierr = DMSetInitialGuess(dms[i],0); 524 if (!mglevels[i+1]->interpolate) { 525 ierr = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr); 526 ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 527 ierr = MatDestroy(&p);CHKERRQ(ierr); 528 } 529 } 530 531 for (i=n-2; i>-1; i--) { 532 ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 533 } 534 ierr = PetscFree(dms);CHKERRQ(ierr); 535 536 /* finest smoother also gets DM but it is not active */ 537 ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 538 ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 539 } 540 541 if (mg->galerkin) { 542 Mat B; 543 mg->galerkinused = PETSC_TRUE; 544 /* currently only handle case where mat and pmat are the same on coarser levels */ 545 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 546 if (!pc->setupcalled) { 547 for (i=n-2; i>-1; i--) { 548 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 549 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 550 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 551 dB = B; 552 } 553 if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 554 } else { 555 for (i=n-2; i>-1; i--) { 556 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 557 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 558 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 559 dB = B; 560 } 561 } 562 } else if (pc->dm && pc->dm->x) { 563 /* need to restrict Jacobian location to coarser meshes for evaluation */ 564 for (i=n-2;i>-1; i--) { 565 if (!mglevels[i]->smoothd->dm->x) { 566 Vec *vecs; 567 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr); 568 mglevels[i]->smoothd->dm->x = vecs[0]; 569 ierr = PetscFree(vecs);CHKERRQ(ierr); 570 } 571 ierr = MatRestrict(mglevels[i+1]->interpolate,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 572 } 573 } 574 575 if (!pc->setupcalled) { 576 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr); 577 578 for (i=0; i<n; i++) { 579 if (monitor) { 580 ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr); 581 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 582 ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 583 } 584 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 585 } 586 for (i=1; i<n; i++) { 587 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 588 if (monitor) { 589 ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr); 590 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 591 ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 592 } 593 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 594 } 595 } 596 for (i=1; i<n; i++) { 597 if (mglevels[i]->restrct && !mglevels[i]->interpolate) { 598 ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr); 599 } 600 if (!mglevels[i]->restrct && mglevels[i]->interpolate) { 601 ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr); 602 } 603 #if defined(PETSC_USE_DEBUG) 604 if (!mglevels[i]->restrct || !mglevels[i]->interpolate) { 605 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 606 } 607 #endif 608 } 609 for (i=0; i<n-1; i++) { 610 if (!mglevels[i]->b) { 611 Vec *vec; 612 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 613 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 614 ierr = VecDestroy(vec);CHKERRQ(ierr); 615 ierr = PetscFree(vec);CHKERRQ(ierr); 616 } 617 if (!mglevels[i]->r && i) { 618 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 619 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 620 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 621 } 622 if (!mglevels[i]->x) { 623 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 624 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 625 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 626 } 627 } 628 if (n != 1 && !mglevels[n-1]->r) { 629 /* PCMGSetR() on the finest level if user did not supply it */ 630 Vec *vec; 631 ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 632 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 633 ierr = VecDestroy(vec);CHKERRQ(ierr); 634 ierr = PetscFree(vec);CHKERRQ(ierr); 635 } 636 } 637 638 639 for (i=1; i<n; i++) { 640 if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 641 /* if doing only down then initial guess is zero */ 642 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 643 } 644 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 645 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 646 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 647 if (!mglevels[i]->residual) { 648 Mat mat; 649 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 650 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 651 } 652 } 653 for (i=1; i<n; i++) { 654 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 655 Mat downmat,downpmat; 656 MatStructure matflag; 657 PetscBool opsset; 658 659 /* check if operators have been set for up, if not use down operators to set them */ 660 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 661 if (!opsset) { 662 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 663 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 664 } 665 666 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 667 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 668 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 669 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 670 } 671 } 672 673 /* 674 If coarse solver is not direct method then DO NOT USE preonly 675 */ 676 ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 677 if (preonly) { 678 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 679 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 680 ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 681 ierr = PetscTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 682 if (!lu && !redundant && !cholesky && !svd) { 683 ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 684 } 685 } 686 687 if (!pc->setupcalled) { 688 if (monitor) { 689 ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr); 690 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 691 ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 692 } 693 ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 694 } 695 696 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 697 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 698 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 699 700 /* 701 Dump the interpolation/restriction matrices plus the 702 Jacobian/stiffness on each level. This allows MATLAB users to 703 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 704 705 Only support one or the other at the same time. 706 */ 707 #if defined(PETSC_USE_SOCKET_VIEWER) 708 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 709 if (dump) { 710 viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 711 } 712 dump = PETSC_FALSE; 713 #endif 714 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 715 if (dump) { 716 viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 717 } 718 719 if (viewer) { 720 for (i=1; i<n; i++) { 721 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 722 } 723 for (i=0; i<n; i++) { 724 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 725 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 726 } 727 } 728 PetscFunctionReturn(0); 729 } 730 731 /* -------------------------------------------------------------------------------------*/ 732 733 #undef __FUNCT__ 734 #define __FUNCT__ "PCMGGetLevels" 735 /*@ 736 PCMGGetLevels - Gets the number of levels to use with MG. 737 738 Not Collective 739 740 Input Parameter: 741 . pc - the preconditioner context 742 743 Output parameter: 744 . levels - the number of levels 745 746 Level: advanced 747 748 .keywords: MG, get, levels, multigrid 749 750 .seealso: PCMGSetLevels() 751 @*/ 752 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 753 { 754 PC_MG *mg = (PC_MG*)pc->data; 755 756 PetscFunctionBegin; 757 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 758 PetscValidIntPointer(levels,2); 759 *levels = mg->nlevels; 760 PetscFunctionReturn(0); 761 } 762 763 #undef __FUNCT__ 764 #define __FUNCT__ "PCMGSetType" 765 /*@ 766 PCMGSetType - Determines the form of multigrid to use: 767 multiplicative, additive, full, or the Kaskade algorithm. 768 769 Logically Collective on PC 770 771 Input Parameters: 772 + pc - the preconditioner context 773 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 774 PC_MG_FULL, PC_MG_KASKADE 775 776 Options Database Key: 777 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 778 additive, full, kaskade 779 780 Level: advanced 781 782 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 783 784 .seealso: PCMGSetLevels() 785 @*/ 786 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 787 { 788 PC_MG *mg = (PC_MG*)pc->data; 789 790 PetscFunctionBegin; 791 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 792 PetscValidLogicalCollectiveEnum(pc,form,2); 793 mg->am = form; 794 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 795 else pc->ops->applyrichardson = 0; 796 PetscFunctionReturn(0); 797 } 798 799 #undef __FUNCT__ 800 #define __FUNCT__ "PCMGSetCycleType" 801 /*@ 802 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 803 complicated cycling. 804 805 Logically Collective on PC 806 807 Input Parameters: 808 + pc - the multigrid context 809 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 810 811 Options Database Key: 812 $ -pc_mg_cycle_type v or w 813 814 Level: advanced 815 816 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 817 818 .seealso: PCMGSetCycleTypeOnLevel() 819 @*/ 820 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 821 { 822 PC_MG *mg = (PC_MG*)pc->data; 823 PC_MG_Levels **mglevels = mg->levels; 824 PetscInt i,levels; 825 826 PetscFunctionBegin; 827 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 828 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 829 PetscValidLogicalCollectiveInt(pc,n,2); 830 levels = mglevels[0]->levels; 831 832 for (i=0; i<levels; i++) { 833 mglevels[i]->cycles = n; 834 } 835 PetscFunctionReturn(0); 836 } 837 838 #undef __FUNCT__ 839 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 840 /*@ 841 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 842 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 843 844 Logically Collective on PC 845 846 Input Parameters: 847 + pc - the multigrid context 848 - n - number of cycles (default is 1) 849 850 Options Database Key: 851 $ -pc_mg_multiplicative_cycles n 852 853 Level: advanced 854 855 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 856 857 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 858 859 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 860 @*/ 861 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 862 { 863 PC_MG *mg = (PC_MG*)pc->data; 864 PC_MG_Levels **mglevels = mg->levels; 865 PetscInt i,levels; 866 867 PetscFunctionBegin; 868 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 869 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 870 PetscValidLogicalCollectiveInt(pc,n,2); 871 levels = mglevels[0]->levels; 872 873 for (i=0; i<levels; i++) { 874 mg->cyclesperpcapply = n; 875 } 876 PetscFunctionReturn(0); 877 } 878 879 #undef __FUNCT__ 880 #define __FUNCT__ "PCMGSetGalerkin" 881 /*@ 882 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 883 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 884 885 Logically Collective on PC 886 887 Input Parameters: 888 . pc - the multigrid context 889 890 Options Database Key: 891 $ -pc_mg_galerkin 892 893 Level: intermediate 894 895 .keywords: MG, set, Galerkin 896 897 .seealso: PCMGGetGalerkin() 898 899 @*/ 900 PetscErrorCode PCMGSetGalerkin(PC pc) 901 { 902 PC_MG *mg = (PC_MG*)pc->data; 903 904 PetscFunctionBegin; 905 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 906 mg->galerkin = PETSC_TRUE; 907 PetscFunctionReturn(0); 908 } 909 910 #undef __FUNCT__ 911 #define __FUNCT__ "PCMGGetGalerkin" 912 /*@ 913 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 914 A_i-1 = r_i * A_i * r_i^t 915 916 Not Collective 917 918 Input Parameter: 919 . pc - the multigrid context 920 921 Output Parameter: 922 . gelerkin - PETSC_TRUE or PETSC_FALSE 923 924 Options Database Key: 925 $ -pc_mg_galerkin 926 927 Level: intermediate 928 929 .keywords: MG, set, Galerkin 930 931 .seealso: PCMGSetGalerkin() 932 933 @*/ 934 PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 935 { 936 PC_MG *mg = (PC_MG*)pc->data; 937 938 PetscFunctionBegin; 939 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 940 *galerkin = mg->galerkin; 941 PetscFunctionReturn(0); 942 } 943 944 #undef __FUNCT__ 945 #define __FUNCT__ "PCMGSetNumberSmoothDown" 946 /*@ 947 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 948 use on all levels. Use PCMGGetSmootherDown() to set different 949 pre-smoothing steps on different levels. 950 951 Logically Collective on PC 952 953 Input Parameters: 954 + mg - the multigrid context 955 - n - the number of smoothing steps 956 957 Options Database Key: 958 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 959 960 Level: advanced 961 962 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 963 964 .seealso: PCMGSetNumberSmoothUp() 965 @*/ 966 PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 967 { 968 PC_MG *mg = (PC_MG*)pc->data; 969 PC_MG_Levels **mglevels = mg->levels; 970 PetscErrorCode ierr; 971 PetscInt i,levels; 972 973 PetscFunctionBegin; 974 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 975 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 976 PetscValidLogicalCollectiveInt(pc,n,2); 977 levels = mglevels[0]->levels; 978 979 for (i=1; i<levels; i++) { 980 /* make sure smoother up and down are different */ 981 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 982 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 983 mg->default_smoothd = n; 984 } 985 PetscFunctionReturn(0); 986 } 987 988 #undef __FUNCT__ 989 #define __FUNCT__ "PCMGSetNumberSmoothUp" 990 /*@ 991 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 992 on all levels. Use PCMGGetSmootherUp() to set different numbers of 993 post-smoothing steps on different levels. 994 995 Logically Collective on PC 996 997 Input Parameters: 998 + mg - the multigrid context 999 - n - the number of smoothing steps 1000 1001 Options Database Key: 1002 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 1003 1004 Level: advanced 1005 1006 Note: this does not set a value on the coarsest grid, since we assume that 1007 there is no separate smooth up on the coarsest grid. 1008 1009 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1010 1011 .seealso: PCMGSetNumberSmoothDown() 1012 @*/ 1013 PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 1014 { 1015 PC_MG *mg = (PC_MG*)pc->data; 1016 PC_MG_Levels **mglevels = mg->levels; 1017 PetscErrorCode ierr; 1018 PetscInt i,levels; 1019 1020 PetscFunctionBegin; 1021 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1022 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1023 PetscValidLogicalCollectiveInt(pc,n,2); 1024 levels = mglevels[0]->levels; 1025 1026 for (i=1; i<levels; i++) { 1027 /* make sure smoother up and down are different */ 1028 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1029 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1030 mg->default_smoothu = n; 1031 } 1032 PetscFunctionReturn(0); 1033 } 1034 1035 /* ----------------------------------------------------------------------------------------*/ 1036 1037 /*MC 1038 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1039 information about the coarser grid matrices and restriction/interpolation operators. 1040 1041 Options Database Keys: 1042 + -pc_mg_levels <nlevels> - number of levels including finest 1043 . -pc_mg_cycles v or w 1044 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 1045 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 1046 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 1047 . -pc_mg_log - log information about time spent on each level of the solver 1048 . -pc_mg_monitor - print information on the multigrid convergence 1049 . -pc_mg_galerkin - use Galerkin process to compute coarser operators 1050 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1051 to the Socket viewer for reading from MATLAB. 1052 1053 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 1054 1055 Level: intermediate 1056 1057 Concepts: multigrid/multilevel 1058 1059 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC 1060 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1061 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1062 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1063 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1064 M*/ 1065 1066 EXTERN_C_BEGIN 1067 #undef __FUNCT__ 1068 #define __FUNCT__ "PCCreate_MG" 1069 PetscErrorCode PCCreate_MG(PC pc) 1070 { 1071 PC_MG *mg; 1072 PetscErrorCode ierr; 1073 1074 PetscFunctionBegin; 1075 ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1076 pc->data = (void*)mg; 1077 mg->nlevels = -1; 1078 1079 pc->ops->apply = PCApply_MG; 1080 pc->ops->setup = PCSetUp_MG; 1081 pc->ops->reset = PCReset_MG; 1082 pc->ops->destroy = PCDestroy_MG; 1083 pc->ops->setfromoptions = PCSetFromOptions_MG; 1084 pc->ops->view = PCView_MG; 1085 PetscFunctionReturn(0); 1086 } 1087 EXTERN_C_END 1088