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