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