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