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 = DMSetFunction(dms[i],0); 487 ierr = DMSetInitialGuess(dms[i],0); 488 if (!mglevels[i+1]->interpolate) { 489 ierr = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr); 490 ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 491 ierr = MatDestroy(&p);CHKERRQ(ierr); 492 } 493 } 494 495 if (!mg->galerkin) { 496 /* each coarse level gets its DM; finest level does not get DM because it shared the outer PC operators */ 497 for (i=n-2; i>-1; i--) { 498 ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 499 } 500 } 501 502 for (i=n-2; i>-1; i--) { 503 ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 504 } 505 ierr = PetscFree(dms);CHKERRQ(ierr); 506 } 507 508 if (mg->galerkin) { 509 Mat B; 510 mg->galerkinused = PETSC_TRUE; 511 /* currently only handle case where mat and pmat are the same on coarser levels */ 512 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 513 if (!pc->setupcalled) { 514 for (i=n-2; i>-1; i--) { 515 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 516 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 517 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 518 dB = B; 519 } 520 if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 521 } else { 522 for (i=n-2; i>-1; i--) { 523 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 524 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 525 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 526 dB = B; 527 } 528 } 529 } 530 531 if (!pc->setupcalled) { 532 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr); 533 534 for (i=0; i<n; i++) { 535 if (monitor) { 536 ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr); 537 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 538 ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 539 } 540 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 541 } 542 for (i=1; i<n; i++) { 543 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 544 if (monitor) { 545 ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr); 546 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 547 ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 548 } 549 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 550 } 551 } 552 for (i=1; i<n; i++) { 553 if (mglevels[i]->restrct && !mglevels[i]->interpolate) { 554 ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr); 555 } 556 if (!mglevels[i]->restrct && mglevels[i]->interpolate) { 557 ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr); 558 } 559 #if defined(PETSC_USE_DEBUG) 560 if (!mglevels[i]->restrct || !mglevels[i]->interpolate) { 561 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 562 } 563 #endif 564 } 565 for (i=0; i<n-1; i++) { 566 if (!mglevels[i]->b) { 567 Vec *vec; 568 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 569 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 570 ierr = VecDestroy(vec);CHKERRQ(ierr); 571 ierr = PetscFree(vec);CHKERRQ(ierr); 572 } 573 if (!mglevels[i]->r && i) { 574 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 575 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 576 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 577 } 578 if (!mglevels[i]->x) { 579 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 580 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 581 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 582 } 583 } 584 if (n != 1 && !mglevels[n-1]->r) { 585 /* PCMGSetR() on the finest level if user did not supply it */ 586 Vec *vec; 587 ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 588 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 589 ierr = VecDestroy(vec);CHKERRQ(ierr); 590 ierr = PetscFree(vec);CHKERRQ(ierr); 591 } 592 } 593 594 595 for (i=1; i<n; i++) { 596 if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 597 /* if doing only down then initial guess is zero */ 598 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 599 } 600 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 601 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 602 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 603 if (!mglevels[i]->residual) { 604 Mat mat; 605 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 606 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 607 } 608 } 609 for (i=1; i<n; i++) { 610 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 611 Mat downmat,downpmat; 612 MatStructure matflag; 613 PetscBool opsset; 614 615 /* check if operators have been set for up, if not use down operators to set them */ 616 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 617 if (!opsset) { 618 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 619 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 620 } 621 622 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 623 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 624 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 625 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 626 } 627 } 628 629 /* 630 If coarse solver is not direct method then DO NOT USE preonly 631 */ 632 ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 633 if (preonly) { 634 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 635 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 636 ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 637 if (!lu && !redundant && !cholesky) { 638 ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 639 } 640 } 641 642 if (!pc->setupcalled) { 643 if (monitor) { 644 ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr); 645 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 646 ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 647 } 648 ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 649 } 650 651 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 652 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 653 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 654 655 /* 656 Dump the interpolation/restriction matrices plus the 657 Jacobian/stiffness on each level. This allows MATLAB users to 658 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 659 660 Only support one or the other at the same time. 661 */ 662 #if defined(PETSC_USE_SOCKET_VIEWER) 663 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 664 if (dump) { 665 viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 666 } 667 dump = PETSC_FALSE; 668 #endif 669 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 670 if (dump) { 671 viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 672 } 673 674 if (viewer) { 675 for (i=1; i<n; i++) { 676 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 677 } 678 for (i=0; i<n; i++) { 679 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 680 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 681 } 682 } 683 PetscFunctionReturn(0); 684 } 685 686 /* -------------------------------------------------------------------------------------*/ 687 688 #undef __FUNCT__ 689 #define __FUNCT__ "PCMGGetLevels" 690 /*@ 691 PCMGGetLevels - Gets the number of levels to use with MG. 692 693 Not Collective 694 695 Input Parameter: 696 . pc - the preconditioner context 697 698 Output parameter: 699 . levels - the number of levels 700 701 Level: advanced 702 703 .keywords: MG, get, levels, multigrid 704 705 .seealso: PCMGSetLevels() 706 @*/ 707 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 708 { 709 PC_MG *mg = (PC_MG*)pc->data; 710 711 PetscFunctionBegin; 712 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 713 PetscValidIntPointer(levels,2); 714 *levels = mg->nlevels; 715 PetscFunctionReturn(0); 716 } 717 718 #undef __FUNCT__ 719 #define __FUNCT__ "PCMGSetType" 720 /*@ 721 PCMGSetType - Determines the form of multigrid to use: 722 multiplicative, additive, full, or the Kaskade algorithm. 723 724 Logically Collective on PC 725 726 Input Parameters: 727 + pc - the preconditioner context 728 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 729 PC_MG_FULL, PC_MG_KASKADE 730 731 Options Database Key: 732 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 733 additive, full, kaskade 734 735 Level: advanced 736 737 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 738 739 .seealso: PCMGSetLevels() 740 @*/ 741 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 742 { 743 PC_MG *mg = (PC_MG*)pc->data; 744 745 PetscFunctionBegin; 746 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 747 PetscValidLogicalCollectiveEnum(pc,form,2); 748 mg->am = form; 749 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 750 else pc->ops->applyrichardson = 0; 751 PetscFunctionReturn(0); 752 } 753 754 #undef __FUNCT__ 755 #define __FUNCT__ "PCMGSetCycleType" 756 /*@ 757 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 758 complicated cycling. 759 760 Logically Collective on PC 761 762 Input Parameters: 763 + pc - the multigrid context 764 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 765 766 Options Database Key: 767 $ -pc_mg_cycle_type v or w 768 769 Level: advanced 770 771 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 772 773 .seealso: PCMGSetCycleTypeOnLevel() 774 @*/ 775 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 776 { 777 PC_MG *mg = (PC_MG*)pc->data; 778 PC_MG_Levels **mglevels = mg->levels; 779 PetscInt i,levels; 780 781 PetscFunctionBegin; 782 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 783 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 784 PetscValidLogicalCollectiveInt(pc,n,2); 785 levels = mglevels[0]->levels; 786 787 for (i=0; i<levels; i++) { 788 mglevels[i]->cycles = n; 789 } 790 PetscFunctionReturn(0); 791 } 792 793 #undef __FUNCT__ 794 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 795 /*@ 796 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 797 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 798 799 Logically Collective on PC 800 801 Input Parameters: 802 + pc - the multigrid context 803 - n - number of cycles (default is 1) 804 805 Options Database Key: 806 $ -pc_mg_multiplicative_cycles n 807 808 Level: advanced 809 810 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 811 812 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 813 814 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 815 @*/ 816 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 817 { 818 PC_MG *mg = (PC_MG*)pc->data; 819 PC_MG_Levels **mglevels = mg->levels; 820 PetscInt i,levels; 821 822 PetscFunctionBegin; 823 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 824 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 825 PetscValidLogicalCollectiveInt(pc,n,2); 826 levels = mglevels[0]->levels; 827 828 for (i=0; i<levels; i++) { 829 mg->cyclesperpcapply = n; 830 } 831 PetscFunctionReturn(0); 832 } 833 834 #undef __FUNCT__ 835 #define __FUNCT__ "PCMGSetGalerkin" 836 /*@ 837 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 838 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 839 840 Logically Collective on PC 841 842 Input Parameters: 843 . pc - the multigrid context 844 845 Options Database Key: 846 $ -pc_mg_galerkin 847 848 Level: intermediate 849 850 .keywords: MG, set, Galerkin 851 852 .seealso: PCMGGetGalerkin() 853 854 @*/ 855 PetscErrorCode PCMGSetGalerkin(PC pc) 856 { 857 PC_MG *mg = (PC_MG*)pc->data; 858 859 PetscFunctionBegin; 860 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 861 mg->galerkin = PETSC_TRUE; 862 PetscFunctionReturn(0); 863 } 864 865 #undef __FUNCT__ 866 #define __FUNCT__ "PCMGGetGalerkin" 867 /*@ 868 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 869 A_i-1 = r_i * A_i * r_i^t 870 871 Not Collective 872 873 Input Parameter: 874 . pc - the multigrid context 875 876 Output Parameter: 877 . gelerkin - PETSC_TRUE or PETSC_FALSE 878 879 Options Database Key: 880 $ -pc_mg_galerkin 881 882 Level: intermediate 883 884 .keywords: MG, set, Galerkin 885 886 .seealso: PCMGSetGalerkin() 887 888 @*/ 889 PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 890 { 891 PC_MG *mg = (PC_MG*)pc->data; 892 893 PetscFunctionBegin; 894 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 895 *galerkin = mg->galerkin; 896 PetscFunctionReturn(0); 897 } 898 899 #undef __FUNCT__ 900 #define __FUNCT__ "PCMGSetNumberSmoothDown" 901 /*@ 902 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 903 use on all levels. Use PCMGGetSmootherDown() to set different 904 pre-smoothing steps on different levels. 905 906 Logically Collective on PC 907 908 Input Parameters: 909 + mg - the multigrid context 910 - n - the number of smoothing steps 911 912 Options Database Key: 913 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 914 915 Level: advanced 916 917 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 918 919 .seealso: PCMGSetNumberSmoothUp() 920 @*/ 921 PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 922 { 923 PC_MG *mg = (PC_MG*)pc->data; 924 PC_MG_Levels **mglevels = mg->levels; 925 PetscErrorCode ierr; 926 PetscInt i,levels; 927 928 PetscFunctionBegin; 929 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 930 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 931 PetscValidLogicalCollectiveInt(pc,n,2); 932 levels = mglevels[0]->levels; 933 934 for (i=1; i<levels; i++) { 935 /* make sure smoother up and down are different */ 936 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 937 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 938 mg->default_smoothd = n; 939 } 940 PetscFunctionReturn(0); 941 } 942 943 #undef __FUNCT__ 944 #define __FUNCT__ "PCMGSetNumberSmoothUp" 945 /*@ 946 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 947 on all levels. Use PCMGGetSmootherUp() to set different numbers of 948 post-smoothing steps on different levels. 949 950 Logically Collective on PC 951 952 Input Parameters: 953 + mg - the multigrid context 954 - n - the number of smoothing steps 955 956 Options Database Key: 957 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 958 959 Level: advanced 960 961 Note: this does not set a value on the coarsest grid, since we assume that 962 there is no separate smooth up on the coarsest grid. 963 964 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 965 966 .seealso: PCMGSetNumberSmoothDown() 967 @*/ 968 PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 969 { 970 PC_MG *mg = (PC_MG*)pc->data; 971 PC_MG_Levels **mglevels = mg->levels; 972 PetscErrorCode ierr; 973 PetscInt i,levels; 974 975 PetscFunctionBegin; 976 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 977 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 978 PetscValidLogicalCollectiveInt(pc,n,2); 979 levels = mglevels[0]->levels; 980 981 for (i=1; i<levels; i++) { 982 /* make sure smoother up and down are different */ 983 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 984 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 985 mg->default_smoothu = n; 986 } 987 PetscFunctionReturn(0); 988 } 989 990 /* ----------------------------------------------------------------------------------------*/ 991 992 /*MC 993 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 994 information about the coarser grid matrices and restriction/interpolation operators. 995 996 Options Database Keys: 997 + -pc_mg_levels <nlevels> - number of levels including finest 998 . -pc_mg_cycles v or w 999 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 1000 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 1001 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 1002 . -pc_mg_log - log information about time spent on each level of the solver 1003 . -pc_mg_monitor - print information on the multigrid convergence 1004 . -pc_mg_galerkin - use Galerkin process to compute coarser operators 1005 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1006 to the Socket viewer for reading from MATLAB. 1007 1008 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 1009 1010 Level: intermediate 1011 1012 Concepts: multigrid/multilevel 1013 1014 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC 1015 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1016 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1017 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1018 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1019 M*/ 1020 1021 EXTERN_C_BEGIN 1022 #undef __FUNCT__ 1023 #define __FUNCT__ "PCCreate_MG" 1024 PetscErrorCode PCCreate_MG(PC pc) 1025 { 1026 PC_MG *mg; 1027 PetscErrorCode ierr; 1028 1029 PetscFunctionBegin; 1030 ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1031 pc->data = (void*)mg; 1032 mg->nlevels = -1; 1033 1034 pc->ops->apply = PCApply_MG; 1035 pc->ops->setup = PCSetUp_MG; 1036 pc->ops->reset = PCReset_MG; 1037 pc->ops->destroy = PCDestroy_MG; 1038 pc->ops->setfromoptions = PCSetFromOptions_MG; 1039 pc->ops->view = PCView_MG; 1040 PetscFunctionReturn(0); 1041 } 1042 EXTERN_C_END 1043