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