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 /* 620 Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level. 621 Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options. 622 But it is safe to use -dm_mat_type. 623 624 The mat type should not be hardcoded like this, we need to find a better way. 625 ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr); 626 */ 627 for (i=n-2; i>-1; i--) { 628 DMKSP kdm; 629 PetscBool dmhasrestrict; 630 ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 631 if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 632 ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 633 /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 634 * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 635 kdm->ops->computerhs = NULL; 636 kdm->rhsctx = NULL; 637 if (!mglevels[i+1]->interpolate) { 638 ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 639 ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 640 if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 641 ierr = VecDestroy(&rscale);CHKERRQ(ierr); 642 ierr = MatDestroy(&p);CHKERRQ(ierr); 643 } 644 ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr); 645 if (dmhasrestrict && !mglevels[i+1]->restrct){ 646 ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr); 647 ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr); 648 ierr = MatDestroy(&p);CHKERRQ(ierr); 649 } 650 } 651 652 for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 653 ierr = PetscFree(dms);CHKERRQ(ierr); 654 } 655 656 if (pc->dm && !pc->setupcalled) { 657 /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 658 ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 659 ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 660 } 661 662 if (mg->galerkin < PC_MG_GALERKIN_NONE) { 663 Mat A,B; 664 PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 665 MatReuse reuse = MAT_INITIAL_MATRIX; 666 667 if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 668 if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE; 669 if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 670 for (i=n-2; i>-1; i--) { 671 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"); 672 if (!mglevels[i+1]->interpolate) { 673 ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 674 } 675 if (!mglevels[i+1]->restrct) { 676 ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 677 } 678 if (reuse == MAT_REUSE_MATRIX) { 679 ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr); 680 } 681 if (doA) { 682 ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr); 683 } 684 if (doB) { 685 ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr); 686 } 687 /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 688 if (!doA && dAeqdB) { 689 if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);} 690 A = B; 691 } else if (!doA && reuse == MAT_INITIAL_MATRIX ) { 692 ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr); 693 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 694 } 695 if (!doB && dAeqdB) { 696 if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);} 697 B = A; 698 } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 699 ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 700 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 701 } 702 if (reuse == MAT_INITIAL_MATRIX) { 703 ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr); 704 ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr); 705 ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr); 706 } 707 dA = A; 708 dB = B; 709 } 710 } 711 if (needRestricts && pc->dm && pc->dm->x) { 712 /* need to restrict Jacobian location to coarser meshes for evaluation */ 713 for (i=n-2; i>-1; i--) { 714 Mat R; 715 Vec rscale; 716 if (!mglevels[i]->smoothd->dm->x) { 717 Vec *vecs; 718 ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 719 mglevels[i]->smoothd->dm->x = vecs[0]; 720 ierr = PetscFree(vecs);CHKERRQ(ierr); 721 } 722 ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 723 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 724 ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 725 ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 726 } 727 } 728 if (needRestricts && pc->dm) { 729 for (i=n-2; i>=0; i--) { 730 DM dmfine,dmcoarse; 731 Mat Restrict,Inject; 732 Vec rscale; 733 ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 734 ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 735 ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 736 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 737 Inject = NULL; /* Callback should create it if it needs Injection */ 738 ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 739 } 740 } 741 742 if (!pc->setupcalled) { 743 for (i=0; i<n; i++) { 744 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 745 } 746 for (i=1; i<n; i++) { 747 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 748 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 749 } 750 } 751 /* insure that if either interpolation or restriction is set the other other one is set */ 752 for (i=1; i<n; i++) { 753 ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr); 754 ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr); 755 } 756 for (i=0; i<n-1; i++) { 757 if (!mglevels[i]->b) { 758 Vec *vec; 759 ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 760 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 761 ierr = VecDestroy(vec);CHKERRQ(ierr); 762 ierr = PetscFree(vec);CHKERRQ(ierr); 763 } 764 if (!mglevels[i]->r && i) { 765 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 766 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 767 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 768 } 769 if (!mglevels[i]->x) { 770 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 771 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 772 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 773 } 774 } 775 if (n != 1 && !mglevels[n-1]->r) { 776 /* PCMGSetR() on the finest level if user did not supply it */ 777 Vec *vec; 778 ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 779 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 780 ierr = VecDestroy(vec);CHKERRQ(ierr); 781 ierr = PetscFree(vec);CHKERRQ(ierr); 782 } 783 } 784 785 if (pc->dm) { 786 /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 787 for (i=0; i<n-1; i++) { 788 if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 789 } 790 } 791 792 for (i=1; i<n; i++) { 793 if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 794 /* if doing only down then initial guess is zero */ 795 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 796 } 797 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 798 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 799 if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 800 pc->failedreason = PC_SUBPC_ERROR; 801 } 802 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 803 if (!mglevels[i]->residual) { 804 Mat mat; 805 ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 806 ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 807 } 808 } 809 for (i=1; i<n; i++) { 810 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 811 Mat downmat,downpmat; 812 813 /* check if operators have been set for up, if not use down operators to set them */ 814 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 815 if (!opsset) { 816 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 817 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 818 } 819 820 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 821 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 822 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 823 if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) { 824 pc->failedreason = PC_SUBPC_ERROR; 825 } 826 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 827 } 828 } 829 830 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 831 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 832 if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 833 pc->failedreason = PC_SUBPC_ERROR; 834 } 835 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 836 837 /* 838 Dump the interpolation/restriction matrices plus the 839 Jacobian/stiffness on each level. This allows MATLAB users to 840 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 841 842 Only support one or the other at the same time. 843 */ 844 #if defined(PETSC_USE_SOCKET_VIEWER) 845 ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 846 if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 847 dump = PETSC_FALSE; 848 #endif 849 ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 850 if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 851 852 if (viewer) { 853 for (i=1; i<n; i++) { 854 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 855 } 856 for (i=0; i<n; i++) { 857 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 858 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 859 } 860 } 861 PetscFunctionReturn(0); 862 } 863 864 /* -------------------------------------------------------------------------------------*/ 865 866 /*@ 867 PCMGGetLevels - Gets the number of levels to use with MG. 868 869 Not Collective 870 871 Input Parameter: 872 . pc - the preconditioner context 873 874 Output parameter: 875 . levels - the number of levels 876 877 Level: advanced 878 879 .keywords: MG, get, levels, multigrid 880 881 .seealso: PCMGSetLevels() 882 @*/ 883 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 884 { 885 PC_MG *mg = (PC_MG*)pc->data; 886 887 PetscFunctionBegin; 888 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 889 PetscValidIntPointer(levels,2); 890 *levels = mg->nlevels; 891 PetscFunctionReturn(0); 892 } 893 894 /*@ 895 PCMGSetType - Determines the form of multigrid to use: 896 multiplicative, additive, full, or the Kaskade algorithm. 897 898 Logically Collective on PC 899 900 Input Parameters: 901 + pc - the preconditioner context 902 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 903 PC_MG_FULL, PC_MG_KASKADE 904 905 Options Database Key: 906 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 907 additive, full, kaskade 908 909 Level: advanced 910 911 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 912 913 .seealso: PCMGSetLevels() 914 @*/ 915 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 916 { 917 PC_MG *mg = (PC_MG*)pc->data; 918 919 PetscFunctionBegin; 920 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 921 PetscValidLogicalCollectiveEnum(pc,form,2); 922 mg->am = form; 923 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 924 else pc->ops->applyrichardson = NULL; 925 PetscFunctionReturn(0); 926 } 927 928 /*@ 929 PCMGGetType - Determines the form of multigrid to use: 930 multiplicative, additive, full, or the Kaskade algorithm. 931 932 Logically Collective on PC 933 934 Input Parameter: 935 . pc - the preconditioner context 936 937 Output Parameter: 938 . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 939 940 941 Level: advanced 942 943 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 944 945 .seealso: PCMGSetLevels() 946 @*/ 947 PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 948 { 949 PC_MG *mg = (PC_MG*)pc->data; 950 951 PetscFunctionBegin; 952 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 953 *type = mg->am; 954 PetscFunctionReturn(0); 955 } 956 957 /*@ 958 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 959 complicated cycling. 960 961 Logically Collective on PC 962 963 Input Parameters: 964 + pc - the multigrid context 965 - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 966 967 Options Database Key: 968 . -pc_mg_cycle_type <v,w> - provide the cycle desired 969 970 Level: advanced 971 972 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 973 974 .seealso: PCMGSetCycleTypeOnLevel() 975 @*/ 976 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 977 { 978 PC_MG *mg = (PC_MG*)pc->data; 979 PC_MG_Levels **mglevels = mg->levels; 980 PetscInt i,levels; 981 982 PetscFunctionBegin; 983 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 984 PetscValidLogicalCollectiveEnum(pc,n,2); 985 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 986 levels = mglevels[0]->levels; 987 for (i=0; i<levels; i++) mglevels[i]->cycles = n; 988 PetscFunctionReturn(0); 989 } 990 991 /*@ 992 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 993 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 994 995 Logically Collective on PC 996 997 Input Parameters: 998 + pc - the multigrid context 999 - n - number of cycles (default is 1) 1000 1001 Options Database Key: 1002 . -pc_mg_multiplicative_cycles n 1003 1004 Level: advanced 1005 1006 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 1007 1008 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 1009 1010 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 1011 @*/ 1012 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 1013 { 1014 PC_MG *mg = (PC_MG*)pc->data; 1015 1016 PetscFunctionBegin; 1017 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1018 PetscValidLogicalCollectiveInt(pc,n,2); 1019 mg->cyclesperpcapply = n; 1020 PetscFunctionReturn(0); 1021 } 1022 1023 PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1024 { 1025 PC_MG *mg = (PC_MG*)pc->data; 1026 1027 PetscFunctionBegin; 1028 mg->galerkin = use; 1029 PetscFunctionReturn(0); 1030 } 1031 1032 /*@ 1033 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 1034 finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1035 1036 Logically Collective on PC 1037 1038 Input Parameters: 1039 + pc - the multigrid context 1040 - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1041 1042 Options Database Key: 1043 . -pc_mg_galerkin <both,pmat,mat,none> 1044 1045 Level: intermediate 1046 1047 Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 1048 use the PCMG construction of the coarser grids. 1049 1050 .keywords: MG, set, Galerkin 1051 1052 .seealso: PCMGGetGalerkin(), PCMGGalerkinType 1053 1054 @*/ 1055 PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1056 { 1057 PetscErrorCode ierr; 1058 1059 PetscFunctionBegin; 1060 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1061 ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr); 1062 PetscFunctionReturn(0); 1063 } 1064 1065 /*@ 1066 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 1067 A_i-1 = r_i * A_i * p_i 1068 1069 Not Collective 1070 1071 Input Parameter: 1072 . pc - the multigrid context 1073 1074 Output Parameter: 1075 . galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL 1076 1077 Level: intermediate 1078 1079 .keywords: MG, set, Galerkin 1080 1081 .seealso: PCMGSetGalerkin(), PCMGGalerkinType 1082 1083 @*/ 1084 PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 1085 { 1086 PC_MG *mg = (PC_MG*)pc->data; 1087 1088 PetscFunctionBegin; 1089 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1090 *galerkin = mg->galerkin; 1091 PetscFunctionReturn(0); 1092 } 1093 1094 /*@ 1095 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 1096 use on all levels. Use PCMGGetSmootherDown() to set different 1097 pre-smoothing steps on different levels. 1098 1099 Logically Collective on PC 1100 1101 Input Parameters: 1102 + mg - the multigrid context 1103 - n - the number of smoothing steps 1104 1105 Options Database Key: 1106 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 1107 1108 Level: advanced 1109 If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the 1110 up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same 1111 number of smoothing steps for pre and post smoothing. 1112 1113 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 1114 1115 .seealso: PCMGSetNumberSmoothUp(), PCMGSetNumberSmooth() 1116 @*/ 1117 PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 1118 { 1119 PC_MG *mg = (PC_MG*)pc->data; 1120 PC_MG_Levels **mglevels = mg->levels; 1121 PetscErrorCode ierr; 1122 PetscInt i,levels; 1123 1124 PetscFunctionBegin; 1125 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1126 PetscValidLogicalCollectiveInt(pc,n,2); 1127 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1128 levels = mglevels[0]->levels; 1129 for (i=1; i<levels; i++) { 1130 PetscInt nc; 1131 ierr = KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);CHKERRQ(ierr); 1132 if (nc == n) continue; 1133 1134 /* make sure smoother up and down are different */ 1135 ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1136 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1137 1138 mg->default_smoothd = n; 1139 } 1140 PetscFunctionReturn(0); 1141 } 1142 1143 /*@ 1144 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 1145 on all levels. Use PCMGGetSmootherUp() to set different numbers of 1146 post-smoothing steps on different levels. 1147 1148 Logically Collective on PC 1149 1150 Input Parameters: 1151 + mg - the multigrid context 1152 - n - the number of smoothing steps 1153 1154 Options Database Key: 1155 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 1156 1157 Level: advanced 1158 1159 Notes: this does not set a value on the coarsest grid, since we assume that 1160 there is no separate smooth up on the coarsest grid. 1161 1162 If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the 1163 up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same 1164 number of smoothing steps for pre and post smoothing. 1165 1166 1167 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1168 1169 .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmooth() 1170 @*/ 1171 PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 1172 { 1173 PC_MG *mg = (PC_MG*)pc->data; 1174 PC_MG_Levels **mglevels = mg->levels; 1175 PetscErrorCode ierr; 1176 PetscInt i,levels; 1177 1178 PetscFunctionBegin; 1179 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1180 PetscValidLogicalCollectiveInt(pc,n,2); 1181 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1182 levels = mglevels[0]->levels; 1183 1184 for (i=1; i<levels; i++) { 1185 if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 1186 PetscInt nc; 1187 ierr = KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);CHKERRQ(ierr); 1188 if (nc == n) continue; 1189 } 1190 1191 /* make sure smoother up and down are different */ 1192 ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1193 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1194 1195 mg->default_smoothu = n; 1196 } 1197 PetscFunctionReturn(0); 1198 } 1199 1200 /*@ 1201 PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1202 on all levels. Use PCMGSetSmoothUp() and PCMGSetSmoothDown() set different numbers of 1203 pre ad post-smoothing steps 1204 1205 Logically Collective on PC 1206 1207 Input Parameters: 1208 + mg - the multigrid context 1209 - n - the number of smoothing steps 1210 1211 Options Database Key: 1212 + -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 1213 . -pc_mg_smooth_down <n> - Sets number of pre-smoothing steps (if setting different pre and post amounts) 1214 - -pc_mg_smooth_up <n> - Sets number of post-smoothing steps (if setting different pre and post amounts) 1215 1216 Level: advanced 1217 1218 Notes: this does not set a value on the coarsest grid, since we assume that 1219 there is no separate smooth up on the coarsest grid. 1220 1221 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1222 1223 .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmoothUp() 1224 @*/ 1225 PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 1226 { 1227 PC_MG *mg = (PC_MG*)pc->data; 1228 PC_MG_Levels **mglevels = mg->levels; 1229 PetscErrorCode ierr; 1230 PetscInt i,levels; 1231 1232 PetscFunctionBegin; 1233 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1234 PetscValidLogicalCollectiveInt(pc,n,2); 1235 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1236 levels = mglevels[0]->levels; 1237 1238 for (i=1; i<levels; i++) { 1239 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1240 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1241 mg->default_smoothu = n; 1242 mg->default_smoothd = n; 1243 } 1244 PetscFunctionReturn(0); 1245 } 1246 1247 /* ----------------------------------------------------------------------------------------*/ 1248 1249 /*MC 1250 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1251 information about the coarser grid matrices and restriction/interpolation operators. 1252 1253 Options Database Keys: 1254 + -pc_mg_levels <nlevels> - number of levels including finest 1255 . -pc_mg_cycle_type <v,w> - 1256 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 1257 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 1258 . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 1259 . -pc_mg_log - log information about time spent on each level of the solver 1260 . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 1261 . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 1262 . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1263 to the Socket viewer for reading from MATLAB. 1264 - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 1265 to the binary output file called binaryoutput 1266 1267 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 1268 1269 When run with a single level the smoother options are used on that level NOT the coarse grid solver options 1270 1271 When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 1272 is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 1273 (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 1274 residual is computed at the end of each cycle. 1275 1276 Level: intermediate 1277 1278 Concepts: multigrid/multilevel 1279 1280 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1281 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1282 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1283 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1284 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1285 M*/ 1286 1287 PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 1288 { 1289 PC_MG *mg; 1290 PetscErrorCode ierr; 1291 1292 PetscFunctionBegin; 1293 ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1294 pc->data = (void*)mg; 1295 mg->nlevels = -1; 1296 mg->am = PC_MG_MULTIPLICATIVE; 1297 mg->galerkin = PC_MG_GALERKIN_NONE; 1298 1299 pc->useAmat = PETSC_TRUE; 1300 1301 pc->ops->apply = PCApply_MG; 1302 pc->ops->setup = PCSetUp_MG; 1303 pc->ops->reset = PCReset_MG; 1304 pc->ops->destroy = PCDestroy_MG; 1305 pc->ops->setfromoptions = PCSetFromOptions_MG; 1306 pc->ops->view = PCView_MG; 1307 1308 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 1309 PetscFunctionReturn(0); 1310 } 1311