1 2 /* 3 Defines the multigrid preconditioner interface. 4 */ 5 #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscpcmg.h" I*/ 6 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "PCMGMCycle_Private" 10 PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 11 { 12 PC_MG *mg = (PC_MG*)pc->data; 13 PC_MG_Levels *mgc,*mglevels = *mglevelsin; 14 PetscErrorCode ierr; 15 PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 16 17 PetscFunctionBegin; 18 if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 19 ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 20 if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 21 if (mglevels->level) { /* not the coarsest grid */ 22 if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 23 ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 24 if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 25 26 /* if on finest level and have convergence criteria set */ 27 if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 28 PetscReal rnorm; 29 ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 30 if (rnorm <= mg->ttol) { 31 if (rnorm < mg->abstol) { 32 *reason = PCRICHARDSON_CONVERGED_ATOL; 33 ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr); 34 } else { 35 *reason = PCRICHARDSON_CONVERGED_RTOL; 36 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); 37 } 38 PetscFunctionReturn(0); 39 } 40 } 41 42 mgc = *(mglevelsin - 1); 43 if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 44 ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 45 if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 46 ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 47 while (cycles--) { 48 ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 49 } 50 if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 51 ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 52 if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 53 if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 54 ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 55 if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 56 } 57 PetscFunctionReturn(0); 58 } 59 60 #undef __FUNCT__ 61 #define __FUNCT__ "PCApplyRichardson_MG" 62 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) 63 { 64 PC_MG *mg = (PC_MG*)pc->data; 65 PC_MG_Levels **mglevels = mg->levels; 66 PetscErrorCode ierr; 67 PetscInt levels = mglevels[0]->levels,i; 68 69 PetscFunctionBegin; 70 mglevels[levels-1]->b = b; 71 mglevels[levels-1]->x = x; 72 73 mg->rtol = rtol; 74 mg->abstol = abstol; 75 mg->dtol = dtol; 76 if (rtol) { 77 /* compute initial residual norm for relative convergence test */ 78 PetscReal rnorm; 79 if (zeroguess) { 80 ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 81 } else { 82 ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 83 ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 84 } 85 mg->ttol = PetscMax(rtol*rnorm,abstol); 86 } else if (abstol) mg->ttol = abstol; 87 else mg->ttol = 0.0; 88 89 /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 90 stop prematurely due to small residual */ 91 for (i=1; i<levels; i++) { 92 ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 93 if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 94 ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 95 } 96 } 97 98 *reason = (PCRichardsonConvergedReason)0; 99 for (i=0; i<its; i++) { 100 ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 101 if (*reason) break; 102 } 103 if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 104 *outits = i; 105 PetscFunctionReturn(0); 106 } 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "PCReset_MG" 110 PetscErrorCode PCReset_MG(PC pc) 111 { 112 PC_MG *mg = (PC_MG*)pc->data; 113 PC_MG_Levels **mglevels = mg->levels; 114 PetscErrorCode ierr; 115 PetscInt i,n; 116 117 PetscFunctionBegin; 118 if (mglevels) { 119 n = mglevels[0]->levels; 120 for (i=0; i<n-1; i++) { 121 ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 122 ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 123 ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 124 ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 125 ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 126 ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 127 } 128 129 for (i=0; i<n; i++) { 130 ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 131 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 132 ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 133 } 134 ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 135 } 136 } 137 PetscFunctionReturn(0); 138 } 139 140 #undef __FUNCT__ 141 #define __FUNCT__ "PCMGSetLevels" 142 /*@C 143 PCMGSetLevels - Sets the number of levels to use with MG. 144 Must be called before any other MG routine. 145 146 Logically Collective on PC 147 148 Input Parameters: 149 + pc - the preconditioner context 150 . levels - the number of levels 151 - comms - optional communicators for each level; this is to allow solving the coarser problems 152 on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 153 154 Level: intermediate 155 156 Notes: 157 If the number of levels is one then the multigrid uses the -mg_levels prefix 158 for setting the level options rather than the -mg_coarse prefix. 159 160 .keywords: MG, set, levels, multigrid 161 162 .seealso: PCMGSetType(), PCMGGetLevels() 163 @*/ 164 PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 165 { 166 PetscErrorCode ierr; 167 PC_MG *mg = (PC_MG*)pc->data; 168 MPI_Comm comm = ((PetscObject)pc)->comm; 169 PC_MG_Levels **mglevels = mg->levels; 170 PetscInt i; 171 PetscMPIInt size; 172 const char *prefix; 173 PC ipc; 174 PetscInt n; 175 176 PetscFunctionBegin; 177 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 178 PetscValidLogicalCollectiveInt(pc,levels,2); 179 if (mg->nlevels == levels) PetscFunctionReturn(0); 180 if (mglevels) { 181 /* changing the number of levels so free up the previous stuff */ 182 ierr = PCReset_MG(pc);CHKERRQ(ierr); 183 n = mglevels[0]->levels; 184 for (i=0; i<n; i++) { 185 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 186 ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 187 } 188 ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 189 ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 190 } 191 ierr = PetscFree(mg->levels);CHKERRQ(ierr); 192 } 193 194 mg->nlevels = levels; 195 196 ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr); 197 ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 198 199 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 200 201 mg->stageApply = 0; 202 for (i=0; i<levels; i++) { 203 ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr); 204 205 mglevels[i]->level = i; 206 mglevels[i]->levels = levels; 207 mglevels[i]->cycles = PC_MG_CYCLE_V; 208 mg->default_smoothu = 2; 209 mg->default_smoothd = 2; 210 mglevels[i]->eventsmoothsetup = 0; 211 mglevels[i]->eventsmoothsolve = 0; 212 mglevels[i]->eventresidual = 0; 213 mglevels[i]->eventinterprestrict = 0; 214 215 if (comms) comm = comms[i]; 216 ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 217 ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 218 ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 219 ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 220 ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 221 ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 222 ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 223 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);CHKERRQ(ierr); 224 ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 225 226 /* do special stuff for coarse grid */ 227 if (!i && levels > 1) { 228 ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 229 230 /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 231 ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 232 ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 233 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 234 if (size > 1) { 235 KSP innerksp; 236 PC innerpc; 237 ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 238 ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr); 239 ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr); 240 ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 241 } else { 242 ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 243 ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 244 } 245 } else { 246 char tprefix[128]; 247 sprintf(tprefix,"mg_levels_%d_",(int)i); 248 ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 249 } 250 ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr); 251 252 mglevels[i]->smoothu = mglevels[i]->smoothd; 253 mg->rtol = 0.0; 254 mg->abstol = 0.0; 255 mg->dtol = 0.0; 256 mg->ttol = 0.0; 257 mg->cyclesperpcapply = 1; 258 } 259 mg->am = PC_MG_MULTIPLICATIVE; 260 mg->levels = mglevels; 261 pc->ops->applyrichardson = PCApplyRichardson_MG; 262 PetscFunctionReturn(0); 263 } 264 265 266 #undef __FUNCT__ 267 #define __FUNCT__ "PCDestroy_MG" 268 PetscErrorCode PCDestroy_MG(PC pc) 269 { 270 PetscErrorCode ierr; 271 PC_MG *mg = (PC_MG*)pc->data; 272 PC_MG_Levels **mglevels = mg->levels; 273 PetscInt i,n; 274 275 PetscFunctionBegin; 276 ierr = PCReset_MG(pc);CHKERRQ(ierr); 277 if (mglevels) { 278 n = mglevels[0]->levels; 279 for (i=0; i<n; i++) { 280 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 281 ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 282 } 283 ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 284 ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 285 } 286 ierr = PetscFree(mg->levels);CHKERRQ(ierr); 287 } 288 ierr = PetscFree(pc->data);CHKERRQ(ierr); 289 PetscFunctionReturn(0); 290 } 291 292 293 294 extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 295 extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 296 extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 297 298 /* 299 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 300 or full cycle of multigrid. 301 302 Note: 303 A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 304 */ 305 #undef __FUNCT__ 306 #define __FUNCT__ "PCApply_MG" 307 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 308 { 309 PC_MG *mg = (PC_MG*)pc->data; 310 PC_MG_Levels **mglevels = mg->levels; 311 PetscErrorCode ierr; 312 PetscInt levels = mglevels[0]->levels,i; 313 314 PetscFunctionBegin; 315 if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 316 /* When the DM is supplying the matrix then it will not exist until here */ 317 for (i=0; i<levels; i++) { 318 if (!mglevels[i]->A) { 319 ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 320 ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 321 } 322 } 323 324 mglevels[levels-1]->b = b; 325 mglevels[levels-1]->x = x; 326 if (mg->am == PC_MG_MULTIPLICATIVE) { 327 ierr = VecSet(x,0.0);CHKERRQ(ierr); 328 for (i=0; i<mg->cyclesperpcapply; i++) { 329 ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr); 330 } 331 } else if (mg->am == PC_MG_ADDITIVE) { 332 ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 333 } else if (mg->am == PC_MG_KASKADE) { 334 ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 335 } else { 336 ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 337 } 338 if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 339 PetscFunctionReturn(0); 340 } 341 342 343 #undef __FUNCT__ 344 #define __FUNCT__ "PCSetFromOptions_MG" 345 PetscErrorCode PCSetFromOptions_MG(PC pc) 346 { 347 PetscErrorCode ierr; 348 PetscInt m,levels = 1,cycles; 349 PetscBool flg,set; 350 PC_MG *mg = (PC_MG*)pc->data; 351 PC_MG_Levels **mglevels = mg->levels; 352 PCMGType mgtype; 353 PCMGCycleType mgctype; 354 355 PetscFunctionBegin; 356 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 357 if (!mg->levels) { 358 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 359 if (!flg && pc->dm) { 360 ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 361 levels++; 362 mg->usedmfornumberoflevels = PETSC_TRUE; 363 } 364 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 365 } 366 mglevels = mg->levels; 367 368 mgctype = (PCMGCycleType) mglevels[0]->cycles; 369 ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 370 if (flg) { 371 ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 372 } 373 flg = PETSC_FALSE; 374 ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr); 375 if (set) { 376 ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr); 377 } 378 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 379 if (flg) { 380 ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 381 } 382 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 383 if (flg) { 384 ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 385 } 386 mgtype = mg->am; 387 ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 388 if (flg) { 389 ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 390 } 391 if (mg->am == PC_MG_MULTIPLICATIVE) { 392 ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 393 if (flg) { 394 ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 395 } 396 } 397 flg = PETSC_FALSE; 398 ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 399 if (flg) { 400 PetscInt i; 401 char eventname[128]; 402 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 403 levels = mglevels[0]->levels; 404 for (i=0; i<levels; i++) { 405 sprintf(eventname,"MGSetup Level %d",(int)i); 406 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 407 sprintf(eventname,"MGSmooth Level %d",(int)i); 408 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 409 if (i) { 410 sprintf(eventname,"MGResid Level %d",(int)i); 411 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 412 sprintf(eventname,"MGInterp Level %d",(int)i); 413 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 414 } 415 } 416 417 #if defined(PETSC_USE_LOG) 418 { 419 const char *sname = "MG Apply"; 420 PetscStageLog stageLog; 421 PetscInt st; 422 423 PetscFunctionBegin; 424 ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 425 for (st = 0; st < stageLog->numStages; ++st) { 426 PetscBool same; 427 428 ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 429 if (same) mg->stageApply = st; 430 } 431 if (!mg->stageApply) { 432 ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 433 } 434 } 435 #endif 436 } 437 ierr = PetscOptionsTail();CHKERRQ(ierr); 438 PetscFunctionReturn(0); 439 } 440 441 const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 442 const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 443 444 #undef __FUNCT__ 445 #define __FUNCT__ "PCView_MG" 446 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 447 { 448 PC_MG *mg = (PC_MG*)pc->data; 449 PC_MG_Levels **mglevels = mg->levels; 450 PetscErrorCode ierr; 451 PetscInt levels = mglevels[0]->levels,i; 452 PetscBool iascii,isbinary,isdraw; 453 454 PetscFunctionBegin; 455 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 456 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 457 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 458 if (iascii) { 459 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); 460 if (mg->am == PC_MG_MULTIPLICATIVE) { 461 ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 462 } 463 if (mg->galerkin) { 464 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 465 } else { 466 ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 467 } 468 for (i=0; i<levels; i++) { 469 if (!i) { 470 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 471 } else { 472 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 473 } 474 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 475 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 476 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 477 if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 478 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 479 } else if (i) { 480 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 481 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 482 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 483 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 484 } 485 } 486 } else if (isbinary) { 487 for (i=levels-1; i>=0; i--) { 488 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 489 if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 490 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 491 } 492 } 493 } else if (isdraw) { 494 PetscDraw draw; 495 PetscReal x,w,y,bottom,th; 496 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 497 ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 498 ierr = PetscDrawStringGetSize(draw,PETSC_NULL,&th);CHKERRQ(ierr); 499 bottom = y - th; 500 for (i=levels-1; i>=0; i--) { 501 if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 502 ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 503 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 504 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 505 } else { 506 w = 0.5*PetscMin(1.0-x,x); 507 ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 508 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 509 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 510 ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 511 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 512 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 513 } 514 ierr = PetscDrawGetBoundingBox(draw,PETSC_NULL,&bottom,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 515 bottom -= th; 516 } 517 } 518 PetscFunctionReturn(0); 519 } 520 521 #include <petsc-private/dmimpl.h> 522 #include <petsc-private/kspimpl.h> 523 524 /* 525 Calls setup for the KSP on each level 526 */ 527 #undef __FUNCT__ 528 #define __FUNCT__ "PCSetUp_MG" 529 PetscErrorCode PCSetUp_MG(PC pc) 530 { 531 PC_MG *mg = (PC_MG*)pc->data; 532 PC_MG_Levels **mglevels = mg->levels; 533 PetscErrorCode ierr; 534 PetscInt i,n = mglevels[0]->levels; 535 PC cpc; 536 PetscBool preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset; 537 Mat dA,dB; 538 MatStructure uflag; 539 Vec tvec; 540 DM *dms; 541 PetscViewer viewer = 0; 542 543 PetscFunctionBegin; 544 /* FIX: Move this to PCSetFromOptions_MG? */ 545 if (mg->usedmfornumberoflevels) { 546 PetscInt levels; 547 ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 548 levels++; 549 if (levels > n) { /* the problem is now being solved on a finer grid */ 550 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 551 n = levels; 552 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 553 mglevels = mg->levels; 554 } 555 } 556 ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 557 558 559 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 560 /* so use those from global PC */ 561 /* Is this what we always want? What if user wants to keep old one? */ 562 ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 563 if (opsset) { 564 Mat mmat; 565 ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr); 566 if (mmat == pc->pmat) opsset = PETSC_FALSE; 567 } 568 if (!opsset) { 569 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); 570 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 571 } 572 573 /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */ 574 if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) { 575 /* construct the interpolation from the DMs */ 576 Mat p; 577 Vec rscale; 578 ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 579 dms[n-1] = pc->dm; 580 for (i=n-2; i>-1; i--) { 581 DMKSP kdm; 582 ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr); 583 ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 584 if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 585 ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 586 /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 587 * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 588 kdm->ops->computerhs = PETSC_NULL; 589 kdm->rhsctx = PETSC_NULL; 590 if (!mglevels[i+1]->interpolate) { 591 ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 592 ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 593 if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 594 ierr = VecDestroy(&rscale);CHKERRQ(ierr); 595 ierr = MatDestroy(&p);CHKERRQ(ierr); 596 } 597 } 598 599 for (i=n-2; i>-1; i--) { 600 ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 601 } 602 ierr = PetscFree(dms);CHKERRQ(ierr); 603 } 604 605 if (pc->dm && !pc->setupcalled) { 606 /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */ 607 ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 608 ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 609 } 610 611 if (mg->galerkin == 1) { 612 Mat B; 613 /* currently only handle case where mat and pmat are the same on coarser levels */ 614 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 615 if (!pc->setupcalled) { 616 for (i=n-2; i>-1; i--) { 617 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 618 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 619 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 620 dB = B; 621 } 622 if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 623 } else { 624 for (i=n-2; i>-1; i--) { 625 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 626 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 627 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 628 dB = B; 629 } 630 } 631 } else if (!mg->galerkin && pc->dm && pc->dm->x) { 632 /* need to restrict Jacobian location to coarser meshes for evaluation */ 633 for (i=n-2; i>-1; i--) { 634 Mat R; 635 Vec rscale; 636 if (!mglevels[i]->smoothd->dm->x) { 637 Vec *vecs; 638 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr); 639 640 mglevels[i]->smoothd->dm->x = vecs[0]; 641 642 ierr = PetscFree(vecs);CHKERRQ(ierr); 643 } 644 ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 645 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 646 ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 647 ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 648 } 649 } 650 if (!mg->galerkin && pc->dm) { 651 for (i=n-2; i>=0; i--) { 652 DM dmfine,dmcoarse; 653 Mat Restrict,Inject; 654 Vec rscale; 655 ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 656 ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 657 ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 658 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 659 Inject = PETSC_NULL; /* Callback should create it if it needs Injection */ 660 ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 661 } 662 } 663 664 if (!pc->setupcalled) { 665 for (i=0; i<n; i++) { 666 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 667 } 668 for (i=1; i<n; i++) { 669 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 670 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 671 } 672 } 673 for (i=1; i<n; i++) { 674 ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 675 ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 676 } 677 for (i=0; i<n-1; i++) { 678 if (!mglevels[i]->b) { 679 Vec *vec; 680 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 681 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 682 ierr = VecDestroy(vec);CHKERRQ(ierr); 683 ierr = PetscFree(vec);CHKERRQ(ierr); 684 } 685 if (!mglevels[i]->r && i) { 686 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 687 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 688 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 689 } 690 if (!mglevels[i]->x) { 691 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 692 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 693 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 694 } 695 } 696 if (n != 1 && !mglevels[n-1]->r) { 697 /* PCMGSetR() on the finest level if user did not supply it */ 698 Vec *vec; 699 ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 700 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 701 ierr = VecDestroy(vec);CHKERRQ(ierr); 702 ierr = PetscFree(vec);CHKERRQ(ierr); 703 } 704 } 705 706 if (pc->dm) { 707 /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 708 for (i=0; i<n-1; i++) { 709 if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 710 } 711 } 712 713 for (i=1; i<n; i++) { 714 if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 715 /* if doing only down then initial guess is zero */ 716 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 717 } 718 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 719 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 720 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 721 if (!mglevels[i]->residual) { 722 Mat mat; 723 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 724 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 725 } 726 } 727 for (i=1; i<n; i++) { 728 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 729 Mat downmat,downpmat; 730 MatStructure matflag; 731 732 /* check if operators have been set for up, if not use down operators to set them */ 733 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 734 if (!opsset) { 735 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 736 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 737 } 738 739 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 740 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 741 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 742 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 743 } 744 } 745 746 /* 747 If coarse solver is not direct method then DO NOT USE preonly 748 */ 749 ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 750 if (preonly) { 751 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 752 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 753 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 754 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 755 if (!lu && !redundant && !cholesky && !svd) { 756 ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 757 } 758 } 759 760 if (!pc->setupcalled) { 761 ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 762 } 763 764 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 765 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 766 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 767 768 /* 769 Dump the interpolation/restriction matrices plus the 770 Jacobian/stiffness on each level. This allows MATLAB users to 771 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 772 773 Only support one or the other at the same time. 774 */ 775 #if defined(PETSC_USE_SOCKET_VIEWER) 776 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 777 if (dump) viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 778 dump = PETSC_FALSE; 779 #endif 780 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 781 if (dump) viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 782 783 if (viewer) { 784 for (i=1; i<n; i++) { 785 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 786 } 787 for (i=0; i<n; i++) { 788 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 789 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 790 } 791 } 792 PetscFunctionReturn(0); 793 } 794 795 /* -------------------------------------------------------------------------------------*/ 796 797 #undef __FUNCT__ 798 #define __FUNCT__ "PCMGGetLevels" 799 /*@ 800 PCMGGetLevels - Gets the number of levels to use with MG. 801 802 Not Collective 803 804 Input Parameter: 805 . pc - the preconditioner context 806 807 Output parameter: 808 . levels - the number of levels 809 810 Level: advanced 811 812 .keywords: MG, get, levels, multigrid 813 814 .seealso: PCMGSetLevels() 815 @*/ 816 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 817 { 818 PC_MG *mg = (PC_MG*)pc->data; 819 820 PetscFunctionBegin; 821 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 822 PetscValidIntPointer(levels,2); 823 *levels = mg->nlevels; 824 PetscFunctionReturn(0); 825 } 826 827 #undef __FUNCT__ 828 #define __FUNCT__ "PCMGSetType" 829 /*@ 830 PCMGSetType - Determines the form of multigrid to use: 831 multiplicative, additive, full, or the Kaskade algorithm. 832 833 Logically Collective on PC 834 835 Input Parameters: 836 + pc - the preconditioner context 837 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 838 PC_MG_FULL, PC_MG_KASKADE 839 840 Options Database Key: 841 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 842 additive, full, kaskade 843 844 Level: advanced 845 846 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 847 848 .seealso: PCMGSetLevels() 849 @*/ 850 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 851 { 852 PC_MG *mg = (PC_MG*)pc->data; 853 854 PetscFunctionBegin; 855 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 856 PetscValidLogicalCollectiveEnum(pc,form,2); 857 mg->am = form; 858 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 859 else pc->ops->applyrichardson = 0; 860 PetscFunctionReturn(0); 861 } 862 863 #undef __FUNCT__ 864 #define __FUNCT__ "PCMGSetCycleType" 865 /*@ 866 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 867 complicated cycling. 868 869 Logically Collective on PC 870 871 Input Parameters: 872 + pc - the multigrid context 873 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 874 875 Options Database Key: 876 $ -pc_mg_cycle_type v or w 877 878 Level: advanced 879 880 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 881 882 .seealso: PCMGSetCycleTypeOnLevel() 883 @*/ 884 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 885 { 886 PC_MG *mg = (PC_MG*)pc->data; 887 PC_MG_Levels **mglevels = mg->levels; 888 PetscInt i,levels; 889 890 PetscFunctionBegin; 891 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 892 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 893 PetscValidLogicalCollectiveInt(pc,n,2); 894 levels = mglevels[0]->levels; 895 896 for (i=0; i<levels; i++) mglevels[i]->cycles = n; 897 PetscFunctionReturn(0); 898 } 899 900 #undef __FUNCT__ 901 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 902 /*@ 903 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 904 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 905 906 Logically Collective on PC 907 908 Input Parameters: 909 + pc - the multigrid context 910 - n - number of cycles (default is 1) 911 912 Options Database Key: 913 $ -pc_mg_multiplicative_cycles n 914 915 Level: advanced 916 917 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 918 919 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 920 921 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 922 @*/ 923 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 924 { 925 PC_MG *mg = (PC_MG*)pc->data; 926 PC_MG_Levels **mglevels = mg->levels; 927 PetscInt i,levels; 928 929 PetscFunctionBegin; 930 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 931 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 932 PetscValidLogicalCollectiveInt(pc,n,2); 933 levels = mglevels[0]->levels; 934 935 for (i=0; i<levels; i++) mg->cyclesperpcapply = n; 936 PetscFunctionReturn(0); 937 } 938 939 #undef __FUNCT__ 940 #define __FUNCT__ "PCMGSetGalerkin" 941 /*@ 942 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 943 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 944 945 Logically Collective on PC 946 947 Input Parameters: 948 + pc - the multigrid context 949 - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 950 951 Options Database Key: 952 $ -pc_mg_galerkin 953 954 Level: intermediate 955 956 .keywords: MG, set, Galerkin 957 958 .seealso: PCMGGetGalerkin() 959 960 @*/ 961 PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 962 { 963 PC_MG *mg = (PC_MG*)pc->data; 964 965 PetscFunctionBegin; 966 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 967 mg->galerkin = use ? 1 : 0; 968 PetscFunctionReturn(0); 969 } 970 971 #undef __FUNCT__ 972 #define __FUNCT__ "PCMGGetGalerkin" 973 /*@ 974 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 975 A_i-1 = r_i * A_i * r_i^t 976 977 Not Collective 978 979 Input Parameter: 980 . pc - the multigrid context 981 982 Output Parameter: 983 . gelerkin - PETSC_TRUE or PETSC_FALSE 984 985 Options Database Key: 986 $ -pc_mg_galerkin 987 988 Level: intermediate 989 990 .keywords: MG, set, Galerkin 991 992 .seealso: PCMGSetGalerkin() 993 994 @*/ 995 PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 996 { 997 PC_MG *mg = (PC_MG*)pc->data; 998 999 PetscFunctionBegin; 1000 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1001 *galerkin = (PetscBool)mg->galerkin; 1002 PetscFunctionReturn(0); 1003 } 1004 1005 #undef __FUNCT__ 1006 #define __FUNCT__ "PCMGSetNumberSmoothDown" 1007 /*@ 1008 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 1009 use on all levels. Use PCMGGetSmootherDown() to set different 1010 pre-smoothing steps on different levels. 1011 1012 Logically Collective on PC 1013 1014 Input Parameters: 1015 + mg - the multigrid context 1016 - n - the number of smoothing steps 1017 1018 Options Database Key: 1019 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 1020 1021 Level: advanced 1022 1023 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 1024 1025 .seealso: PCMGSetNumberSmoothUp() 1026 @*/ 1027 PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 1028 { 1029 PC_MG *mg = (PC_MG*)pc->data; 1030 PC_MG_Levels **mglevels = mg->levels; 1031 PetscErrorCode ierr; 1032 PetscInt i,levels; 1033 1034 PetscFunctionBegin; 1035 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1036 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1037 PetscValidLogicalCollectiveInt(pc,n,2); 1038 levels = mglevels[0]->levels; 1039 1040 for (i=1; i<levels; i++) { 1041 /* make sure smoother up and down are different */ 1042 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1043 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1044 1045 mg->default_smoothd = n; 1046 } 1047 PetscFunctionReturn(0); 1048 } 1049 1050 #undef __FUNCT__ 1051 #define __FUNCT__ "PCMGSetNumberSmoothUp" 1052 /*@ 1053 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 1054 on all levels. Use PCMGGetSmootherUp() to set different numbers of 1055 post-smoothing steps on different levels. 1056 1057 Logically Collective on PC 1058 1059 Input Parameters: 1060 + mg - the multigrid context 1061 - n - the number of smoothing steps 1062 1063 Options Database Key: 1064 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 1065 1066 Level: advanced 1067 1068 Note: this does not set a value on the coarsest grid, since we assume that 1069 there is no separate smooth up on the coarsest grid. 1070 1071 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1072 1073 .seealso: PCMGSetNumberSmoothDown() 1074 @*/ 1075 PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 1076 { 1077 PC_MG *mg = (PC_MG*)pc->data; 1078 PC_MG_Levels **mglevels = mg->levels; 1079 PetscErrorCode ierr; 1080 PetscInt i,levels; 1081 1082 PetscFunctionBegin; 1083 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1084 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1085 PetscValidLogicalCollectiveInt(pc,n,2); 1086 levels = mglevels[0]->levels; 1087 1088 for (i=1; i<levels; i++) { 1089 /* make sure smoother up and down are different */ 1090 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1091 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1092 1093 mg->default_smoothu = n; 1094 } 1095 PetscFunctionReturn(0); 1096 } 1097 1098 /* ----------------------------------------------------------------------------------------*/ 1099 1100 /*MC 1101 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1102 information about the coarser grid matrices and restriction/interpolation operators. 1103 1104 Options Database Keys: 1105 + -pc_mg_levels <nlevels> - number of levels including finest 1106 . -pc_mg_cycles v or w 1107 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 1108 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 1109 . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 1110 . -pc_mg_log - log information about time spent on each level of the solver 1111 . -pc_mg_monitor - print information on the multigrid convergence 1112 . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 1113 . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 1114 . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1115 to the Socket viewer for reading from MATLAB. 1116 - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 1117 to the binary output file called binaryoutput 1118 1119 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 1120 1121 When run with a single level the smoother options are used on that level NOT the coarse grid solver options 1122 1123 Level: intermediate 1124 1125 Concepts: multigrid/multilevel 1126 1127 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1128 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1129 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1130 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1131 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1132 M*/ 1133 1134 EXTERN_C_BEGIN 1135 #undef __FUNCT__ 1136 #define __FUNCT__ "PCCreate_MG" 1137 PetscErrorCode PCCreate_MG(PC pc) 1138 { 1139 PC_MG *mg; 1140 PetscErrorCode ierr; 1141 1142 PetscFunctionBegin; 1143 ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1144 pc->data = (void*)mg; 1145 mg->nlevels = -1; 1146 1147 pc->ops->apply = PCApply_MG; 1148 pc->ops->setup = PCSetUp_MG; 1149 pc->ops->reset = PCReset_MG; 1150 pc->ops->destroy = PCDestroy_MG; 1151 pc->ops->setfromoptions = PCSetFromOptions_MG; 1152 pc->ops->view = PCView_MG; 1153 PetscFunctionReturn(0); 1154 } 1155 EXTERN_C_END 1156