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