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