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