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