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