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