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