1 2 /* 3 Defines the multigrid preconditioner interface. 4 */ 5 #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscksp.h" I*/ 6 #include <petscdm.h> 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 #include <petscdraw.h> 446 #undef __FUNCT__ 447 #define __FUNCT__ "PCView_MG" 448 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 449 { 450 PC_MG *mg = (PC_MG*)pc->data; 451 PC_MG_Levels **mglevels = mg->levels; 452 PetscErrorCode ierr; 453 PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 454 PetscBool iascii,isbinary,isdraw; 455 456 PetscFunctionBegin; 457 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 458 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 459 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 460 if (iascii) { 461 const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 462 ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 463 if (mg->am == PC_MG_MULTIPLICATIVE) { 464 ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 465 } 466 if (mg->galerkin) { 467 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 468 } else { 469 ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 470 } 471 for (i=0; i<levels; i++) { 472 if (!i) { 473 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 474 } else { 475 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 476 } 477 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 478 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 479 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 480 if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 481 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 482 } else if (i) { 483 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 484 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 485 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 486 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 487 } 488 } 489 } else if (isbinary) { 490 for (i=levels-1; i>=0; i--) { 491 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 492 if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 493 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 494 } 495 } 496 } else if (isdraw) { 497 PetscDraw draw; 498 PetscReal x,w,y,bottom,th; 499 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 500 ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 501 ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 502 bottom = y - th; 503 for (i=levels-1; i>=0; i--) { 504 if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 505 ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 506 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 507 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 508 } else { 509 w = 0.5*PetscMin(1.0-x,x); 510 ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 511 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 512 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 513 ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 514 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 515 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 516 } 517 ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 518 bottom -= th; 519 } 520 } 521 PetscFunctionReturn(0); 522 } 523 524 #include <petsc-private/dmimpl.h> 525 #include <petsc-private/kspimpl.h> 526 527 /* 528 Calls setup for the KSP on each level 529 */ 530 #undef __FUNCT__ 531 #define __FUNCT__ "PCSetUp_MG" 532 PetscErrorCode PCSetUp_MG(PC pc) 533 { 534 PC_MG *mg = (PC_MG*)pc->data; 535 PC_MG_Levels **mglevels = mg->levels; 536 PetscErrorCode ierr; 537 PetscInt i,n = mglevels[0]->levels; 538 PC cpc; 539 PetscBool preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset,use_amat; 540 Mat dA,dB; 541 MatStructure uflag; 542 Vec tvec; 543 DM *dms; 544 PetscViewer viewer = 0; 545 546 PetscFunctionBegin; 547 /* FIX: Move this to PCSetFromOptions_MG? */ 548 if (mg->usedmfornumberoflevels) { 549 PetscInt levels; 550 ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 551 levels++; 552 if (levels > n) { /* the problem is now being solved on a finer grid */ 553 ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 554 n = levels; 555 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 556 mglevels = mg->levels; 557 } 558 } 559 ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 560 561 562 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 563 /* so use those from global PC */ 564 /* Is this what we always want? What if user wants to keep old one? */ 565 ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 566 if (opsset) { 567 Mat mmat; 568 ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat,NULL);CHKERRQ(ierr); 569 if (mmat == pc->pmat) opsset = PETSC_FALSE; 570 } 571 572 if (!opsset) { 573 ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 574 if(use_amat){ 575 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); 576 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 577 } 578 else { 579 ierr = PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr); 580 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat,pc->flag);CHKERRQ(ierr); 581 } 582 } 583 584 /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */ 585 if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) { 586 /* construct the interpolation from the DMs */ 587 Mat p; 588 Vec rscale; 589 ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 590 dms[n-1] = pc->dm; 591 for (i=n-2; i>-1; i--) { 592 DMKSP kdm; 593 ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr); 594 ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 595 if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 596 ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 597 /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 598 * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 599 kdm->ops->computerhs = NULL; 600 kdm->rhsctx = NULL; 601 if (!mglevels[i+1]->interpolate) { 602 ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 603 ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 604 if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 605 ierr = VecDestroy(&rscale);CHKERRQ(ierr); 606 ierr = MatDestroy(&p);CHKERRQ(ierr); 607 } 608 } 609 610 for (i=n-2; i>-1; i--) { 611 ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 612 } 613 ierr = PetscFree(dms);CHKERRQ(ierr); 614 } 615 616 if (pc->dm && !pc->setupcalled) { 617 /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */ 618 ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 619 ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 620 } 621 622 if (mg->galerkin == 1) { 623 Mat B; 624 /* currently only handle case where mat and pmat are the same on coarser levels */ 625 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 626 if (!pc->setupcalled) { 627 for (i=n-2; i>-1; i--) { 628 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 629 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 630 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 631 dB = B; 632 } 633 if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 634 } else { 635 for (i=n-2; i>-1; i--) { 636 ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B,NULL);CHKERRQ(ierr); 637 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 638 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 639 dB = B; 640 } 641 } 642 } else if (!mg->galerkin && pc->dm && pc->dm->x) { 643 /* need to restrict Jacobian location to coarser meshes for evaluation */ 644 for (i=n-2; i>-1; i--) { 645 Mat R; 646 Vec rscale; 647 if (!mglevels[i]->smoothd->dm->x) { 648 Vec *vecs; 649 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 650 651 mglevels[i]->smoothd->dm->x = vecs[0]; 652 653 ierr = PetscFree(vecs);CHKERRQ(ierr); 654 } 655 ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 656 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 657 ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 658 ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 659 } 660 } 661 if (!mg->galerkin && pc->dm) { 662 for (i=n-2; i>=0; i--) { 663 DM dmfine,dmcoarse; 664 Mat Restrict,Inject; 665 Vec rscale; 666 ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 667 ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 668 ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 669 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 670 Inject = NULL; /* Callback should create it if it needs Injection */ 671 ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 672 } 673 } 674 675 if (!pc->setupcalled) { 676 for (i=0; i<n; i++) { 677 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 678 } 679 for (i=1; i<n; i++) { 680 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 681 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 682 } 683 } 684 for (i=1; i<n; i++) { 685 ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 686 ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 687 } 688 for (i=0; i<n-1; i++) { 689 if (!mglevels[i]->b) { 690 Vec *vec; 691 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 692 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 693 ierr = VecDestroy(vec);CHKERRQ(ierr); 694 ierr = PetscFree(vec);CHKERRQ(ierr); 695 } 696 if (!mglevels[i]->r && i) { 697 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 698 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 699 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 700 } 701 if (!mglevels[i]->x) { 702 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 703 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 704 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 705 } 706 } 707 if (n != 1 && !mglevels[n-1]->r) { 708 /* PCMGSetR() on the finest level if user did not supply it */ 709 Vec *vec; 710 ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 711 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 712 ierr = VecDestroy(vec);CHKERRQ(ierr); 713 ierr = PetscFree(vec);CHKERRQ(ierr); 714 } 715 } 716 717 if (pc->dm) { 718 /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 719 for (i=0; i<n-1; i++) { 720 if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 721 } 722 } 723 724 for (i=1; i<n; i++) { 725 if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 726 /* if doing only down then initial guess is zero */ 727 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 728 } 729 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 730 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 731 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 732 if (!mglevels[i]->residual) { 733 Mat mat; 734 ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat,NULL);CHKERRQ(ierr); 735 ierr = PCMGSetResidual(pc,i,PCMGResidual_Default,mat);CHKERRQ(ierr); 736 } 737 } 738 for (i=1; i<n; i++) { 739 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 740 Mat downmat,downpmat; 741 MatStructure matflag; 742 743 /* check if operators have been set for up, if not use down operators to set them */ 744 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 745 if (!opsset) { 746 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 747 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 748 } 749 750 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 751 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 752 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 753 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 754 } 755 } 756 757 /* 758 If coarse solver is not direct method then DO NOT USE preonly 759 */ 760 ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 761 if (preonly) { 762 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 763 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 764 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 765 ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 766 if (!lu && !redundant && !cholesky && !svd) { 767 ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 768 } 769 } 770 771 if (!pc->setupcalled) { 772 ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 773 } 774 775 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 776 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 777 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 778 779 /* 780 Dump the interpolation/restriction matrices plus the 781 Jacobian/stiffness on each level. This allows MATLAB users to 782 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 783 784 Only support one or the other at the same time. 785 */ 786 #if defined(PETSC_USE_SOCKET_VIEWER) 787 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 788 if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 789 dump = PETSC_FALSE; 790 #endif 791 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 792 if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 793 794 if (viewer) { 795 for (i=1; i<n; i++) { 796 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 797 } 798 for (i=0; i<n; i++) { 799 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 800 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 801 } 802 } 803 PetscFunctionReturn(0); 804 } 805 806 /* -------------------------------------------------------------------------------------*/ 807 808 #undef __FUNCT__ 809 #define __FUNCT__ "PCMGGetLevels" 810 /*@ 811 PCMGGetLevels - Gets the number of levels to use with MG. 812 813 Not Collective 814 815 Input Parameter: 816 . pc - the preconditioner context 817 818 Output parameter: 819 . levels - the number of levels 820 821 Level: advanced 822 823 .keywords: MG, get, levels, multigrid 824 825 .seealso: PCMGSetLevels() 826 @*/ 827 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 828 { 829 PC_MG *mg = (PC_MG*)pc->data; 830 831 PetscFunctionBegin; 832 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 833 PetscValidIntPointer(levels,2); 834 *levels = mg->nlevels; 835 PetscFunctionReturn(0); 836 } 837 838 #undef __FUNCT__ 839 #define __FUNCT__ "PCMGSetType" 840 /*@ 841 PCMGSetType - Determines the form of multigrid to use: 842 multiplicative, additive, full, or the Kaskade algorithm. 843 844 Logically Collective on PC 845 846 Input Parameters: 847 + pc - the preconditioner context 848 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 849 PC_MG_FULL, PC_MG_KASKADE 850 851 Options Database Key: 852 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 853 additive, full, kaskade 854 855 Level: advanced 856 857 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 858 859 .seealso: PCMGSetLevels() 860 @*/ 861 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 862 { 863 PC_MG *mg = (PC_MG*)pc->data; 864 865 PetscFunctionBegin; 866 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 867 PetscValidLogicalCollectiveEnum(pc,form,2); 868 mg->am = form; 869 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 870 else pc->ops->applyrichardson = 0; 871 PetscFunctionReturn(0); 872 } 873 874 #undef __FUNCT__ 875 #define __FUNCT__ "PCMGSetCycleType" 876 /*@ 877 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 878 complicated cycling. 879 880 Logically Collective on PC 881 882 Input Parameters: 883 + pc - the multigrid context 884 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 885 886 Options Database Key: 887 $ -pc_mg_cycle_type v or w 888 889 Level: advanced 890 891 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 892 893 .seealso: PCMGSetCycleTypeOnLevel() 894 @*/ 895 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 896 { 897 PC_MG *mg = (PC_MG*)pc->data; 898 PC_MG_Levels **mglevels = mg->levels; 899 PetscInt i,levels; 900 901 PetscFunctionBegin; 902 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 903 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 904 PetscValidLogicalCollectiveInt(pc,n,2); 905 levels = mglevels[0]->levels; 906 907 for (i=0; i<levels; i++) mglevels[i]->cycles = n; 908 PetscFunctionReturn(0); 909 } 910 911 #undef __FUNCT__ 912 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 913 /*@ 914 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 915 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 916 917 Logically Collective on PC 918 919 Input Parameters: 920 + pc - the multigrid context 921 - n - number of cycles (default is 1) 922 923 Options Database Key: 924 $ -pc_mg_multiplicative_cycles n 925 926 Level: advanced 927 928 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 929 930 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 931 932 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 933 @*/ 934 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 935 { 936 PC_MG *mg = (PC_MG*)pc->data; 937 PC_MG_Levels **mglevels = mg->levels; 938 PetscInt i,levels; 939 940 PetscFunctionBegin; 941 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 942 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 943 PetscValidLogicalCollectiveInt(pc,n,2); 944 levels = mglevels[0]->levels; 945 946 for (i=0; i<levels; i++) mg->cyclesperpcapply = n; 947 PetscFunctionReturn(0); 948 } 949 950 #undef __FUNCT__ 951 #define __FUNCT__ "PCMGSetGalerkin" 952 /*@ 953 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 954 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 955 956 Logically Collective on PC 957 958 Input Parameters: 959 + pc - the multigrid context 960 - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 961 962 Options Database Key: 963 $ -pc_mg_galerkin 964 965 Level: intermediate 966 967 .keywords: MG, set, Galerkin 968 969 .seealso: PCMGGetGalerkin() 970 971 @*/ 972 PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 973 { 974 PC_MG *mg = (PC_MG*)pc->data; 975 976 PetscFunctionBegin; 977 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 978 mg->galerkin = use ? 1 : 0; 979 PetscFunctionReturn(0); 980 } 981 982 #undef __FUNCT__ 983 #define __FUNCT__ "PCMGGetGalerkin" 984 /*@ 985 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 986 A_i-1 = r_i * A_i * r_i^t 987 988 Not Collective 989 990 Input Parameter: 991 . pc - the multigrid context 992 993 Output Parameter: 994 . gelerkin - PETSC_TRUE or PETSC_FALSE 995 996 Options Database Key: 997 $ -pc_mg_galerkin 998 999 Level: intermediate 1000 1001 .keywords: MG, set, Galerkin 1002 1003 .seealso: PCMGSetGalerkin() 1004 1005 @*/ 1006 PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 1007 { 1008 PC_MG *mg = (PC_MG*)pc->data; 1009 1010 PetscFunctionBegin; 1011 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1012 *galerkin = (PetscBool)mg->galerkin; 1013 PetscFunctionReturn(0); 1014 } 1015 1016 #undef __FUNCT__ 1017 #define __FUNCT__ "PCMGSetNumberSmoothDown" 1018 /*@ 1019 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 1020 use on all levels. Use PCMGGetSmootherDown() to set different 1021 pre-smoothing steps on different levels. 1022 1023 Logically Collective on PC 1024 1025 Input Parameters: 1026 + mg - the multigrid context 1027 - n - the number of smoothing steps 1028 1029 Options Database Key: 1030 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 1031 1032 Level: advanced 1033 1034 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 1035 1036 .seealso: PCMGSetNumberSmoothUp() 1037 @*/ 1038 PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 1039 { 1040 PC_MG *mg = (PC_MG*)pc->data; 1041 PC_MG_Levels **mglevels = mg->levels; 1042 PetscErrorCode ierr; 1043 PetscInt i,levels; 1044 1045 PetscFunctionBegin; 1046 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1047 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1048 PetscValidLogicalCollectiveInt(pc,n,2); 1049 levels = mglevels[0]->levels; 1050 1051 for (i=1; i<levels; i++) { 1052 /* make sure smoother up and down are different */ 1053 ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1054 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1055 1056 mg->default_smoothd = n; 1057 } 1058 PetscFunctionReturn(0); 1059 } 1060 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "PCMGSetNumberSmoothUp" 1063 /*@ 1064 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 1065 on all levels. Use PCMGGetSmootherUp() to set different numbers of 1066 post-smoothing steps on different levels. 1067 1068 Logically Collective on PC 1069 1070 Input Parameters: 1071 + mg - the multigrid context 1072 - n - the number of smoothing steps 1073 1074 Options Database Key: 1075 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 1076 1077 Level: advanced 1078 1079 Note: this does not set a value on the coarsest grid, since we assume that 1080 there is no separate smooth up on the coarsest grid. 1081 1082 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1083 1084 .seealso: PCMGSetNumberSmoothDown() 1085 @*/ 1086 PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 1087 { 1088 PC_MG *mg = (PC_MG*)pc->data; 1089 PC_MG_Levels **mglevels = mg->levels; 1090 PetscErrorCode ierr; 1091 PetscInt i,levels; 1092 1093 PetscFunctionBegin; 1094 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1095 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1096 PetscValidLogicalCollectiveInt(pc,n,2); 1097 levels = mglevels[0]->levels; 1098 1099 for (i=1; i<levels; i++) { 1100 /* make sure smoother up and down are different */ 1101 ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1102 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1103 1104 mg->default_smoothu = n; 1105 } 1106 PetscFunctionReturn(0); 1107 } 1108 1109 /* ----------------------------------------------------------------------------------------*/ 1110 1111 /*MC 1112 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1113 information about the coarser grid matrices and restriction/interpolation operators. 1114 1115 Options Database Keys: 1116 + -pc_mg_levels <nlevels> - number of levels including finest 1117 . -pc_mg_cycles v or w 1118 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 1119 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 1120 . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 1121 . -pc_mg_log - log information about time spent on each level of the solver 1122 . -pc_mg_monitor - print information on the multigrid convergence 1123 . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 1124 . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 1125 . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1126 to the Socket viewer for reading from MATLAB. 1127 - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 1128 to the binary output file called binaryoutput 1129 1130 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 1131 1132 When run with a single level the smoother options are used on that level NOT the coarse grid solver options 1133 1134 Level: intermediate 1135 1136 Concepts: multigrid/multilevel 1137 1138 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1139 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1140 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1141 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1142 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1143 M*/ 1144 1145 #undef __FUNCT__ 1146 #define __FUNCT__ "PCCreate_MG" 1147 PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 1148 { 1149 PC_MG *mg; 1150 PetscErrorCode ierr; 1151 1152 PetscFunctionBegin; 1153 ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1154 pc->data = (void*)mg; 1155 mg->nlevels = -1; 1156 1157 pc->useAmat = PETSC_TRUE; 1158 1159 pc->ops->apply = PCApply_MG; 1160 pc->ops->setup = PCSetUp_MG; 1161 pc->ops->reset = PCReset_MG; 1162 pc->ops->destroy = PCDestroy_MG; 1163 pc->ops->setfromoptions = PCSetFromOptions_MG; 1164 pc->ops->view = PCView_MG; 1165 PetscFunctionReturn(0); 1166 } 1167