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