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