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